xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision fb8e56e08d4d0bfe9fc63603ca1f7fddd68abbdb)
1 #include <petsc-private/linesearchimpl.h> /*I "petscsnes.h" I*/
2 
3 PetscBool         SNESLineSearchRegisterAllCalled = PETSC_FALSE;
4 PetscFunctionList SNESLineSearchList              = NULL;
5 
6 PetscClassId  SNESLINESEARCH_CLASSID;
7 PetscLogEvent SNESLineSearch_Apply;
8 
9 #undef __FUNCT__
10 #define __FUNCT__ "SNESLineSearchCreate"
11 /*@
12    SNESLineSearchCreate - Creates the line search context.
13 
14    Logically Collective on Comm
15 
16    Input Parameters:
17 .  comm - MPI communicator for the line search (typically from the associated SNES context).
18 
19    Output Parameters:
20 .  outlinesearch - the new linesearch context
21 
22    Level: beginner
23 
24 .keywords: LineSearch, create, context
25 
26 .seealso: LineSearchDestroy()
27 @*/
28 
29 PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch)
30 {
31   PetscErrorCode ierr;
32   SNESLineSearch linesearch;
33 
34   PetscFunctionBegin;
35   PetscValidPointer(outlinesearch,2);
36   *outlinesearch = NULL;
37 
38   ierr = PetscHeaderCreate(linesearch,_p_LineSearch,struct _LineSearchOps,SNESLINESEARCH_CLASSID, "SNESLineSearch","Linesearch","SNESLineSearch",comm,SNESLineSearchDestroy,SNESLineSearchView);CHKERRQ(ierr);
39 
40   linesearch->ops->precheck  = NULL;
41   linesearch->ops->postcheck = NULL;
42 
43   linesearch->vec_sol_new  = NULL;
44   linesearch->vec_func_new = NULL;
45   linesearch->vec_sol      = NULL;
46   linesearch->vec_func     = NULL;
47   linesearch->vec_update   = NULL;
48 
49   linesearch->lambda       = 1.0;
50   linesearch->fnorm        = 1.0;
51   linesearch->ynorm        = 1.0;
52   linesearch->xnorm        = 1.0;
53   linesearch->success      = PETSC_TRUE;
54   linesearch->norms        = PETSC_TRUE;
55   linesearch->keeplambda   = PETSC_FALSE;
56   linesearch->damping      = 1.0;
57   linesearch->maxstep      = 1e8;
58   linesearch->steptol      = 1e-12;
59   linesearch->rtol         = 1e-8;
60   linesearch->atol         = 1e-15;
61   linesearch->ltol         = 1e-8;
62   linesearch->precheckctx  = NULL;
63   linesearch->postcheckctx = NULL;
64   linesearch->max_its      = 1;
65   linesearch->setupcalled  = PETSC_FALSE;
66   *outlinesearch           = linesearch;
67   PetscFunctionReturn(0);
68 }
69 
70 #undef __FUNCT__
71 #define __FUNCT__ "SNESLineSearchSetUp"
72 /*@
73    SNESLineSearchSetUp - Prepares the line search for being applied by allocating
74    any required vectors.
75 
76    Collective on SNESLineSearch
77 
78    Input Parameters:
79 .  linesearch - The LineSearch instance.
80 
81    Notes:
82    For most cases, this needn't be called outside of SNESLineSearchApply().
83    The only current case where this is called outside of this is for the VI
84    solvers, which modify the solution and work vectors before the first call
85    of SNESLineSearchApply, requiring the SNESLineSearch work vectors to be
86    allocated upfront.
87 
88    Level: advanced
89 
90 .keywords: SNESLineSearch, SetUp
91 
92 .seealso: SNESLineSearchReset()
93 @*/
94 
95 PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch)
96 {
97   PetscErrorCode ierr;
98 
99   PetscFunctionBegin;
100   if (!((PetscObject)linesearch)->type_name) {
101     ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr);
102   }
103   if (!linesearch->setupcalled) {
104     if (!linesearch->vec_sol_new) {
105       ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new);CHKERRQ(ierr);
106     }
107     if (!linesearch->vec_func_new) {
108       ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new);CHKERRQ(ierr);
109     }
110     if (linesearch->ops->setup) {
111       ierr = (*linesearch->ops->setup)(linesearch);CHKERRQ(ierr);
112     }
113     if (!linesearch->ops->snesfunc) {ierr = SNESLineSearchSetFunction(linesearch,SNESComputeFunction);CHKERRQ(ierr);}
114     linesearch->lambda      = linesearch->damping;
115     linesearch->setupcalled = PETSC_TRUE;
116   }
117   PetscFunctionReturn(0);
118 }
119 
120 #undef __FUNCT__
121 #define __FUNCT__ "SNESLineSearchReset"
122 
123 /*@
124    SNESLineSearchReset - Undoes the SetUp and deletes any Vecs or Mats allocated by the line search.
125 
126    Collective on SNESLineSearch
127 
128    Input Parameters:
129 .  linesearch - The LineSearch instance.
130 
131    Level: intermediate
132 
133 .keywords: SNESLineSearch, Reset
134 
135 .seealso: SNESLineSearchSetUp()
136 @*/
137 
138 PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch)
139 {
140   PetscErrorCode ierr;
141 
142   PetscFunctionBegin;
143   if (linesearch->ops->reset) (*linesearch->ops->reset)(linesearch);
144 
145   ierr = VecDestroy(&linesearch->vec_sol_new);CHKERRQ(ierr);
146   ierr = VecDestroy(&linesearch->vec_func_new);CHKERRQ(ierr);
147 
148   ierr = VecDestroyVecs(linesearch->nwork, &linesearch->work);CHKERRQ(ierr);
149 
150   linesearch->nwork       = 0;
151   linesearch->setupcalled = PETSC_FALSE;
152   PetscFunctionReturn(0);
153 }
154 
155 #undef __FUNCT__
156 #define __FUNCT__ "SNESLineSearchSetFunction"
157 /*@C
158    SNESLineSearchSetPreCheck - Sets the function evaluation context on the
159 
160    Input Parameters:
161 .  linesearch - the SNESLineSearch context
162 
163    Output Parameters:
164 +  func       - [optional] function evaluation routine
165 
166    Level: developer
167 
168 .keywords: get, linesearch, pre-check
169 
170 .seealso: SNESLineSearchGetPostCheck(), SNESLineSearchSetPreCheck()
171 @*/
172 PetscErrorCode  SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES,Vec,Vec))
173 {
174   PetscFunctionBegin;
175   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
176   linesearch->ops->snesfunc = func;
177   PetscFunctionReturn(0);
178 }
179 
180 
181 /*MC
182     SNESLineSearchPreCheckFunction - functional form passed to check before line search is called
183 
184      Synopsis:
185      #include "petscsnes.h"
186      SNESLineSearchPreCheckFunction(SNESLineSearch snes,Vec x,Vec y, PetscBool *changed);
187 
188        Input Parameters:
189 +      x - solution vector
190 .      y - search direction vector
191 -      changed - flag to indicate the precheck changed x or y.
192 
193    Level: advanced
194 
195 .seealso:   SNESLineSearchSetPreCheck(), SNESLineSearchSetPostCheck()
196 M*/
197 
198 #undef __FUNCT__
199 #define __FUNCT__ "SNESLineSearchSetPreCheck"
200 /*@C
201    SNESLineSearchSetPreCheck - Sets a pre-check function for the line search routine.
202 
203    Logically Collective on SNESLineSearch
204 
205    Input Parameters:
206 +  linesearch - the SNESLineSearch context
207 .  SNESLineSearchPreCheckFunction - [optional] function evaluation routine
208 -  ctx        - [optional] user-defined context for private data for the
209                 function evaluation routine (may be NULL)
210 
211 
212    Level: intermediate
213 
214 .keywords: set, linesearch, pre-check
215 
216 .seealso: SNESLineSearchSetPostCheck()
217 @*/
218 PetscErrorCode  SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*SNESLineSearchPreCheckFunction)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void *ctx)
219 {
220   PetscFunctionBegin;
221   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
222   if (SNESLineSearchPreCheckFunction) linesearch->ops->precheck = SNESLineSearchPreCheckFunction;
223   if (ctx) linesearch->precheckctx = ctx;
224   PetscFunctionReturn(0);
225 }
226 
227 #undef __FUNCT__
228 #define __FUNCT__ "SNESLineSearchGetPreCheck"
229 /*@C
230    SNESLineSearchSetPreCheck - Gets the pre-check function for the line search routine.
231 
232    Input Parameters:
233 .  linesearch - the SNESLineSearch context
234 
235    Output Parameters:
236 +  func       - [optional] function evaluation routine
237 -  ctx        - [optional] user-defined context for private data for the
238                 function evaluation routine (may be NULL)
239 
240    Level: intermediate
241 
242 .keywords: get, linesearch, pre-check
243 
244 .seealso: SNESLineSearchGetPostCheck(), SNESLineSearchSetPreCheck()
245 @*/
246 PetscErrorCode  SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**SNESLineSearchPreCheckFunction)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void **ctx)
247 {
248   PetscFunctionBegin;
249   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
250   if (SNESLineSearchPreCheckFunction) *SNESLineSearchPreCheckFunction = linesearch->ops->precheck;
251   if (ctx) *ctx = linesearch->precheckctx;
252   PetscFunctionReturn(0);
253 }
254 
255 /*MC
256     SNESLineSearchPostheckFunction - functional form that is called after line search is complete
257 
258      Synopsis:
259      #include "petscsnes.h"
260      SNESLineSearchPostheckFunction(SNESLineSearch linesearch,Vec x,Vec y,  Vec w, *changed_y, PetscBool *changed_w);
261 
262      Input Parameters:
263 +      x - old solution vector
264 .      y - search direction vector
265 .      w - new solution vector
266 .      changed_y - indicates that the line search changed y
267 -      changed_w - indicates that the line search changed w
268 
269    Level: advanced
270 
271 .seealso:   SNESLineSearchSetPreCheck(), SNESLineSearchSetPostCheck()
272 M*/
273 
274 #undef __FUNCT__
275 #define __FUNCT__ "SNESLineSearchSetPostCheck"
276 /*@C
277    SNESLineSearchSetPostCheck - Sets a post-check function for the line search routine.
278 
279    Logically Collective on SNESLineSearch
280 
281    Input Parameters:
282 +  linesearch - the SNESLineSearch context
283 .  SNESLineSearchPostCheckFunction - [optional] function evaluation routine
284 -  ctx        - [optional] user-defined context for private data for the
285                 function evaluation routine (may be NULL)
286 
287    Level: intermediate
288 
289 .keywords: set, linesearch, post-check
290 
291 .seealso: SNESLineSearchSetPreCheck()
292 @*/
293 PetscErrorCode  SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*SNESLineSearchPostCheckFunction)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void *ctx)
294 {
295   PetscFunctionBegin;
296   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
297   if (SNESLineSearchPostCheckFunction) linesearch->ops->postcheck = SNESLineSearchPostCheckFunction;
298   if (ctx) linesearch->postcheckctx = ctx;
299   PetscFunctionReturn(0);
300 }
301 
302 #undef __FUNCT__
303 #define __FUNCT__ "SNESLineSearchGetPostCheck"
304 /*@C
305    SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine.
306 
307    Input Parameters:
308 .  linesearch - the SNESLineSearch context
309 
310    Output Parameters:
311 +  SNESLineSearchPostCheckFunction - [optional] function evaluation routine
312 -  ctx        - [optional] user-defined context for private data for the
313                 function evaluation routine (may be NULL)
314 
315    Level: intermediate
316 
317 .keywords: get, linesearch, post-check
318 
319 .seealso: SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck()
320 @*/
321 PetscErrorCode  SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**SNESLineSearchPostCheckFunction)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void **ctx)
322 {
323   PetscFunctionBegin;
324   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
325   if (SNESLineSearchPostCheckFunction) *SNESLineSearchPostCheckFunction = linesearch->ops->postcheck;
326   if (ctx) *ctx = linesearch->postcheckctx;
327   PetscFunctionReturn(0);
328 }
329 
330 #undef __FUNCT__
331 #define __FUNCT__ "SNESLineSearchPreCheck"
332 /*@
333    SNESLineSearchPreCheck - Prepares the line search for being applied.
334 
335    Logically Collective on SNESLineSearch
336 
337    Input Parameters:
338 +  linesearch - The linesearch instance.
339 .  X - The current solution
340 -  Y - The step direction
341 
342    Output Parameters:
343 .  changed - Indicator that the precheck routine has changed anything
344 
345    Level: Beginner
346 
347 .keywords: SNESLineSearch, Create
348 
349 .seealso: SNESLineSearchPostCheck()
350 @*/
351 PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed)
352 {
353   PetscErrorCode ierr;
354 
355   PetscFunctionBegin;
356   *changed = PETSC_FALSE;
357   if (linesearch->ops->precheck) {
358     ierr = (*linesearch->ops->precheck)(linesearch, X, Y, changed, linesearch->precheckctx);CHKERRQ(ierr);
359     PetscValidLogicalCollectiveBool(linesearch,*changed,4);
360   }
361   PetscFunctionReturn(0);
362 }
363 
364 #undef __FUNCT__
365 #define __FUNCT__ "SNESLineSearchPostCheck"
366 /*@
367    SNESLineSearchPostCheck - Prepares the line search for being applied.
368 
369    Logically Collective on SNESLineSearch
370 
371    Input Parameters:
372 +  linesearch - The linesearch context
373 .  X - The last solution
374 .  Y - The step direction
375 -  W - The updated solution, W = X + lambda*Y for some lambda
376 
377    Output Parameters:
378 +  changed_Y - Indicator if the direction Y has been changed.
379 -  changed_W - Indicator if the new candidate solution W has been changed.
380 
381    Level: Intermediate
382 
383 .keywords: SNESLineSearch, Create
384 
385 .seealso: SNESLineSearchPreCheck()
386 @*/
387 PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch,Vec X,Vec Y,Vec W,PetscBool *changed_Y,PetscBool *changed_W)
388 {
389   PetscErrorCode ierr;
390 
391   PetscFunctionBegin;
392   *changed_Y = PETSC_FALSE;
393   *changed_W = PETSC_FALSE;
394   if (linesearch->ops->postcheck) {
395     ierr = (*linesearch->ops->postcheck)(linesearch,X,Y,W,changed_Y,changed_W,linesearch->postcheckctx);CHKERRQ(ierr);
396     PetscValidLogicalCollectiveBool(linesearch,*changed_Y,5);
397     PetscValidLogicalCollectiveBool(linesearch,*changed_W,6);
398   }
399   PetscFunctionReturn(0);
400 }
401 
402 #undef __FUNCT__
403 #define __FUNCT__ "SNESLineSearchPreCheckPicard"
404 /*@C
405    SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration
406 
407    Logically Collective on SNESLineSearch
408 
409    Input Arguments:
410 +  linesearch - linesearch context
411 .  X - base state for this step
412 .  Y - initial correction
413 
414    Output Arguments:
415 +  Y - correction, possibly modified
416 -  changed - flag indicating that Y was modified
417 
418    Options Database Key:
419 +  -snes_linesearch_precheck_picard - activate this routine
420 -  -snes_linesearch_precheck_picard_angle - angle
421 
422    Level: advanced
423 
424    Notes:
425    This function should be passed to SNESLineSearchSetPreCheck()
426 
427    The justification for this method involves the linear convergence of a Picard iteration
428    so the Picard linearization should be provided in place of the "Jacobian". This correction
429    is generally not useful when using a Newton linearization.
430 
431    Reference:
432    Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology.
433 
434 .seealso: SNESLineSearchSetPreCheck()
435 @*/
436 PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed,void *ctx)
437 {
438   PetscErrorCode ierr;
439   PetscReal      angle = *(PetscReal*)linesearch->precheckctx;
440   Vec            Ylast;
441   PetscScalar    dot;
442   PetscInt       iter;
443   PetscReal      ynorm,ylastnorm,theta,angle_radians;
444   SNES           snes;
445 
446   PetscFunctionBegin;
447   ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
448   ierr = PetscObjectQuery((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject*)&Ylast);CHKERRQ(ierr);
449   if (!Ylast) {
450     ierr = VecDuplicate(Y,&Ylast);CHKERRQ(ierr);
451     ierr = PetscObjectCompose((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject)Ylast);CHKERRQ(ierr);
452     ierr = PetscObjectDereference((PetscObject)Ylast);CHKERRQ(ierr);
453   }
454   ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr);
455   if (iter < 2) {
456     ierr     = VecCopy(Y,Ylast);CHKERRQ(ierr);
457     *changed = PETSC_FALSE;
458     PetscFunctionReturn(0);
459   }
460 
461   ierr = VecDot(Y,Ylast,&dot);CHKERRQ(ierr);
462   ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr);
463   ierr = VecNorm(Ylast,NORM_2,&ylastnorm);CHKERRQ(ierr);
464   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
465   theta         = acos((double)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0));
466   angle_radians = angle * PETSC_PI / 180.;
467   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
468     /* Modify the step Y */
469     PetscReal alpha,ydiffnorm;
470     ierr  = VecAXPY(Ylast,-1.0,Y);CHKERRQ(ierr);
471     ierr  = VecNorm(Ylast,NORM_2,&ydiffnorm);CHKERRQ(ierr);
472     alpha = ylastnorm / ydiffnorm;
473     ierr  = VecCopy(Y,Ylast);CHKERRQ(ierr);
474     ierr  = VecScale(Y,alpha);CHKERRQ(ierr);
475     ierr  = PetscInfo3(snes,"Angle %14.12e degrees less than threshold %14.12e, corrected step by alpha=%14.12e\n",(double)(theta*180./PETSC_PI),(double)angle,(double)alpha);CHKERRQ(ierr);
476   } else {
477     ierr     = PetscInfo2(snes,"Angle %14.12e degrees exceeds threshold %14.12e, no correction applied\n",(double)(theta*180./PETSC_PI),(double)angle);CHKERRQ(ierr);
478     ierr     = VecCopy(Y,Ylast);CHKERRQ(ierr);
479     *changed = PETSC_FALSE;
480   }
481   PetscFunctionReturn(0);
482 }
483 
484 #undef __FUNCT__
485 #define __FUNCT__ "SNESLineSearchApply"
486 /*@
487    SNESLineSearchApply - Computes the line-search update.
488 
489    Collective on SNESLineSearch
490 
491    Input Parameters:
492 +  linesearch - The linesearch context
493 .  X - The current solution
494 .  F - The current function
495 .  fnorm - The current norm
496 -  Y - The search direction
497 
498    Output Parameters:
499 +  X - The new solution
500 .  F - The new function
501 -  fnorm - The new function norm
502 
503    Options Database Keys:
504 + -snes_linesearch_type - basic, bt, l2, cp, shell
505 . -snes_linesearch_monitor - Print progress of line searches
506 . -snes_linesearch_damping - The linesearch damping parameter
507 . -snes_linesearch_norms   - Turn on/off the linesearch norms
508 . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess
509 - -snes_linesearch_max_it - The number of iterations for iterative line searches
510 
511    Notes:
512    This is typically called from within a SNESSolve() implementation in order to
513    help with convergence of the nonlinear method.  Various SNES types use line searches
514    in different ways, but the overarching theme is that a line search is used to determine
515    an optimal damping parameter of a step at each iteration of the method.  Each
516    application of the line search may invoke SNESComputeFunction several times, and
517    therefore may be fairly expensive.
518 
519    Level: Intermediate
520 
521 .keywords: SNESLineSearch, Create
522 
523 .seealso: SNESLineSearchCreate(), SNESLineSearchPreCheck(), SNESLineSearchPostCheck(), SNESSolve(), SNESComputeFunction()
524 @*/
525 PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y)
526 {
527   PetscErrorCode ierr;
528 
529   PetscFunctionBegin;
530   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
531   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
532   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
533   PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
534 
535   linesearch->success = PETSC_TRUE;
536 
537   linesearch->vec_sol    = X;
538   linesearch->vec_update = Y;
539   linesearch->vec_func   = F;
540 
541   ierr = SNESLineSearchSetUp(linesearch);CHKERRQ(ierr);
542 
543   if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */
544 
545   if (fnorm) linesearch->fnorm = *fnorm;
546   else {
547     ierr = VecNorm(F, NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
548   }
549 
550   ierr = PetscLogEventBegin(SNESLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
551 
552   ierr = (*linesearch->ops->apply)(linesearch);CHKERRQ(ierr);
553 
554   ierr = PetscLogEventEnd(SNESLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
555 
556   if (fnorm) *fnorm = linesearch->fnorm;
557   PetscFunctionReturn(0);
558 }
559 
560 #undef __FUNCT__
561 #define __FUNCT__ "SNESLineSearchDestroy"
562 /*@
563    SNESLineSearchDestroy - Destroys the line search instance.
564 
565    Collective on SNESLineSearch
566 
567    Input Parameters:
568 .  linesearch - The linesearch context
569 
570    Level: Intermediate
571 
572 .keywords: SNESLineSearch, Destroy
573 
574 .seealso: SNESLineSearchCreate(), SNESLineSearchReset(), SNESDestroy()
575 @*/
576 PetscErrorCode SNESLineSearchDestroy(SNESLineSearch * linesearch)
577 {
578   PetscErrorCode ierr;
579 
580   PetscFunctionBegin;
581   if (!*linesearch) PetscFunctionReturn(0);
582   PetscValidHeaderSpecific((*linesearch),SNESLINESEARCH_CLASSID,1);
583   if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = 0; PetscFunctionReturn(0);}
584   ierr = PetscObjectAMSViewOff((PetscObject)*linesearch);CHKERRQ(ierr);
585   ierr = SNESLineSearchReset(*linesearch);CHKERRQ(ierr);
586   if ((*linesearch)->ops->destroy) (*linesearch)->ops->destroy(*linesearch);
587   ierr = PetscViewerDestroy(&(*linesearch)->monitor);CHKERRQ(ierr);
588   ierr = PetscHeaderDestroy(linesearch);CHKERRQ(ierr);
589   PetscFunctionReturn(0);
590 }
591 
592 #undef __FUNCT__
593 #define __FUNCT__ "SNESLineSearchSetMonitor"
594 /*@
595    SNESLineSearchSetMonitor - Turns on/off printing useful information and debugging output about the line search.
596 
597    Input Parameters:
598 +  snes - nonlinear context obtained from SNESCreate()
599 -  flg - PETSC_TRUE to monitor the line search
600 
601    Logically Collective on SNES
602 
603    Options Database:
604 .   -snes_linesearch_monitor - enables the monitor
605 
606    Level: intermediate
607 
608 
609 .seealso: SNESLineSearchGetMonitor(), PetscViewer
610 @*/
611 PetscErrorCode  SNESLineSearchSetMonitor(SNESLineSearch linesearch, PetscBool flg)
612 {
613   PetscErrorCode ierr;
614 
615   PetscFunctionBegin;
616   if (flg && !linesearch->monitor) {
617     ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)linesearch),"stdout",&linesearch->monitor);CHKERRQ(ierr);
618   } else if (!flg && linesearch->monitor) {
619     ierr = PetscViewerDestroy(&linesearch->monitor);CHKERRQ(ierr);
620   }
621   PetscFunctionReturn(0);
622 }
623 
624 #undef __FUNCT__
625 #define __FUNCT__ "SNESLineSearchGetMonitor"
626 /*@
627    SNESLineSearchGetMonitor - Gets the PetscViewer instance for the line search monitor.
628 
629    Input Parameters:
630 .  linesearch - linesearch context
631 
632    Input Parameters:
633 .  monitor - monitor context
634 
635    Logically Collective on SNES
636 
637 
638    Options Database Keys:
639 .   -snes_linesearch_monitor - enables the monitor
640 
641    Level: intermediate
642 
643 
644 .seealso: SNESLineSearchSetMonitor(), PetscViewer
645 @*/
646 PetscErrorCode  SNESLineSearchGetMonitor(SNESLineSearch linesearch, PetscViewer *monitor)
647 {
648   PetscFunctionBegin;
649   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
650   if (monitor) {
651     PetscValidPointer(monitor, 2);
652     *monitor = linesearch->monitor;
653   }
654   PetscFunctionReturn(0);
655 }
656 
657 #undef __FUNCT__
658 #define __FUNCT__ "SNESLineSearchSetFromOptions"
659 /*@
660    SNESLineSearchSetFromOptions - Sets options for the line search
661 
662    Input Parameters:
663 .  linesearch - linesearch context
664 
665    Options Database Keys:
666 + -snes_linesearch_type <type> - basic, bt, l2, cp, shell
667 . -snes_linesearch_order <order> - 1, 2, 3.  Most types only support certain orders (bt supports 2 or 3)
668 . -snes_linesearch_norms   - Turn on/off the linesearch norms for the basic linesearch type
669 . -snes_linesearch_minlambda - The minimum step length
670 . -snes_linesearch_maxstep - The maximum step size
671 . -snes_linesearch_rtol - Relative tolerance for iterative line searches
672 . -snes_linesearch_atol - Absolute tolerance for iterative line searches
673 . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches
674 . -snes_linesearch_max_it - The number of iterations for iterative line searches
675 . -snes_linesearch_monitor - Print progress of line searches
676 . -snes_linesearch_damping - The linesearch damping parameter
677 . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess.
678 . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method
679 - -snes_linesearch_precheck_picard_angle - Angle used in picard precheck method
680 
681    Logically Collective on SNESLineSearch
682 
683    Level: intermediate
684 
685 .seealso: SNESLineSearchCreate(), SNESLineSearchSetOrder(), SNESLineSearchSetType(), SNESLineSearchSetTolerances(), SNESLineSearchSetDamping(), SNESLineSearchPreCheckPicard()
686 @*/
687 PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch)
688 {
689   PetscErrorCode ierr;
690   const char     *deft = SNESLINESEARCHBASIC;
691   char           type[256];
692   PetscBool      flg, set;
693 
694   PetscFunctionBegin;
695   if (!SNESLineSearchRegisterAllCalled) {ierr = SNESLineSearchRegisterAll();CHKERRQ(ierr);}
696 
697   ierr = PetscObjectOptionsBegin((PetscObject)linesearch);CHKERRQ(ierr);
698   if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name;
699   ierr = PetscOptionsList("-snes_linesearch_type","Linesearch type","SNESLineSearchSetType",SNESLineSearchList,deft,type,256,&flg);CHKERRQ(ierr);
700   if (flg) {
701     ierr = SNESLineSearchSetType(linesearch,type);CHKERRQ(ierr);
702   } else if (!((PetscObject)linesearch)->type_name) {
703     ierr = SNESLineSearchSetType(linesearch,deft);CHKERRQ(ierr);
704   }
705 
706   ierr = PetscOptionsBool("-snes_linesearch_monitor","Print progress of line searches","SNESSNESLineSearchSetMonitor",
707                           linesearch->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
708   if (set) {ierr = SNESLineSearchSetMonitor(linesearch,flg);CHKERRQ(ierr);}
709 
710   /* tolerances */
711   ierr = PetscOptionsReal("-snes_linesearch_minlambda","Minimum step length","SNESLineSearchSetTolerances",linesearch->steptol,&linesearch->steptol,0);CHKERRQ(ierr);
712   ierr = PetscOptionsReal("-snes_linesearch_maxstep","Maximum step size","SNESLineSearchSetTolerances",linesearch->maxstep,&linesearch->maxstep,0);CHKERRQ(ierr);
713   ierr = PetscOptionsReal("-snes_linesearch_rtol","Relative tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->rtol,&linesearch->rtol,0);CHKERRQ(ierr);
714   ierr = PetscOptionsReal("-snes_linesearch_atol","Absolute tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->atol,&linesearch->atol,0);CHKERRQ(ierr);
715   ierr = PetscOptionsReal("-snes_linesearch_ltol","Change in lambda tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->ltol,&linesearch->ltol,0);CHKERRQ(ierr);
716   ierr = PetscOptionsInt("-snes_linesearch_max_it","Maximum iterations for iterative line searches","SNESLineSearchSetTolerances",linesearch->max_its,&linesearch->max_its,0);CHKERRQ(ierr);
717 
718   /* damping parameters */
719   ierr = PetscOptionsReal("-snes_linesearch_damping","Line search damping and initial step guess","SNESLineSearchSetDamping",linesearch->damping,&linesearch->damping,0);CHKERRQ(ierr);
720 
721   ierr = PetscOptionsBool("-snes_linesearch_keeplambda","Use previous lambda as damping","SNESLineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,0);CHKERRQ(ierr);
722 
723   /* precheck */
724   ierr = PetscOptionsBool("-snes_linesearch_precheck_picard","Use a correction that sometimes improves convergence of Picard iteration","SNESLineSearchPreCheckPicard",flg,&flg,&set);CHKERRQ(ierr);
725   if (set) {
726     if (flg) {
727       linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */
728 
729       ierr = PetscOptionsReal("-snes_linesearch_precheck_picard_angle","Maximum angle at which to activate the correction",
730                               "none",linesearch->precheck_picard_angle,&linesearch->precheck_picard_angle,NULL);CHKERRQ(ierr);
731       ierr = SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&linesearch->precheck_picard_angle);CHKERRQ(ierr);
732     } else {
733       ierr = SNESLineSearchSetPreCheck(linesearch,NULL,NULL);CHKERRQ(ierr);
734     }
735   }
736   ierr = PetscOptionsInt("-snes_linesearch_order","Order of approximation used in the line search","SNESLineSearchSetOrder",linesearch->order,&linesearch->order,0);CHKERRQ(ierr);
737   ierr = PetscOptionsBool("-snes_linesearch_norms","Compute final norms in line search","SNESLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,0);CHKERRQ(ierr);
738 
739   if (linesearch->ops->setfromoptions) {
740     (*linesearch->ops->setfromoptions)(linesearch);CHKERRQ(ierr);
741   }
742 
743   ierr = PetscObjectProcessOptionsHandlers((PetscObject)linesearch);CHKERRQ(ierr);
744   ierr = PetscOptionsEnd();CHKERRQ(ierr);
745   PetscFunctionReturn(0);
746 }
747 
748 #undef __FUNCT__
749 #define __FUNCT__ "SNESLineSearchView"
750 /*@
751    SNESLineSearchView - Prints useful information about the line search not
752    related to an individual call.
753 
754    Input Parameters:
755 .  linesearch - linesearch context
756 
757    Logically Collective on SNESLineSearch
758 
759    Level: intermediate
760 
761 .seealso: SNESLineSearchCreate()
762 @*/
763 PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer)
764 {
765   PetscErrorCode ierr;
766   PetscBool      iascii;
767 
768   PetscFunctionBegin;
769   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
770   if (!viewer) {
771     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch),&viewer);CHKERRQ(ierr);
772   }
773   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
774   PetscCheckSameComm(linesearch,1,viewer,2);
775 
776   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
777   if (iascii) {
778     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)linesearch,viewer,"SNESLineSearch Object");CHKERRQ(ierr);
779     if (linesearch->ops->view) {
780       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
781       ierr = (*linesearch->ops->view)(linesearch,viewer);CHKERRQ(ierr);
782       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
783     }
784     ierr = PetscViewerASCIIPrintf(viewer,"  maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep,(double)linesearch->steptol);CHKERRQ(ierr);
785     ierr = PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol,(double)linesearch->atol,(double)linesearch->ltol);CHKERRQ(ierr);
786     ierr = PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D\n", linesearch->max_its);CHKERRQ(ierr);
787     if (linesearch->ops->precheck) {
788       if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) {
789         ierr = PetscViewerASCIIPrintf(viewer,"  using precheck step to speed up Picard convergence\n", linesearch->max_its);CHKERRQ(ierr);
790       } else {
791         ierr = PetscViewerASCIIPrintf(viewer,"  using user-defined precheck step\n", linesearch->max_its);CHKERRQ(ierr);
792       }
793     }
794     if (linesearch->ops->postcheck) {
795       ierr = PetscViewerASCIIPrintf(viewer,"  using user-defined postcheck step\n", linesearch->max_its);CHKERRQ(ierr);
796     }
797   }
798   PetscFunctionReturn(0);
799 }
800 
801 #undef __FUNCT__
802 #define __FUNCT__ "SNESLineSearchSetType"
803 /*@C
804    SNESLineSearchSetType - Sets the linesearch type
805 
806    Input Parameters:
807 +  linesearch - linesearch context
808 -  type - The type of line search to be used
809 
810    Available Types:
811 +  basic - Simple damping line search.
812 .  bt - Backtracking line search over the L2 norm of the function
813 .  l2 - Secant line search over the L2 norm of the function
814 .  cp - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x)
815 -  shell - User provided SNESLineSearch implementation
816 
817    Logically Collective on SNESLineSearch
818 
819    Level: intermediate
820 
821 
822 .seealso: SNESLineSearchCreate()
823 @*/
824 PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type)
825 {
826   PetscErrorCode ierr,(*r)(SNESLineSearch);
827   PetscBool      match;
828 
829   PetscFunctionBegin;
830   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
831   PetscValidCharPointer(type,2);
832 
833   ierr = PetscObjectTypeCompare((PetscObject)linesearch,type,&match);CHKERRQ(ierr);
834   if (match) PetscFunctionReturn(0);
835 
836   ierr = PetscFunctionListFind(SNESLineSearchList,type,&r);CHKERRQ(ierr);
837   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type);
838   /* Destroy the previous private linesearch context */
839   if (linesearch->ops->destroy) {
840     ierr = (*(linesearch)->ops->destroy)(linesearch);CHKERRQ(ierr);
841 
842     linesearch->ops->destroy = NULL;
843   }
844   /* Reinitialize function pointers in SNESLineSearchOps structure */
845   linesearch->ops->apply          = 0;
846   linesearch->ops->view           = 0;
847   linesearch->ops->setfromoptions = 0;
848   linesearch->ops->destroy        = 0;
849 
850   ierr = PetscObjectChangeTypeName((PetscObject)linesearch,type);CHKERRQ(ierr);
851   ierr = (*r)(linesearch);CHKERRQ(ierr);
852   PetscFunctionReturn(0);
853 }
854 
855 #undef __FUNCT__
856 #define __FUNCT__ "SNESLineSearchSetSNES"
857 /*@
858    SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation.
859 
860    Input Parameters:
861 +  linesearch - linesearch context
862 -  snes - The snes instance
863 
864    Level: developer
865 
866    Notes:
867    This happens automatically when the line search is gotten/created with
868    SNESGetLineSearch().  This routine is therefore mainly called within SNES
869    implementations.
870 
871    Level: developer
872 
873 .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
874 @*/
875 PetscErrorCode  SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes)
876 {
877   PetscFunctionBegin;
878   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
879   PetscValidHeaderSpecific(snes,SNES_CLASSID,2);
880   linesearch->snes = snes;
881   PetscFunctionReturn(0);
882 }
883 
884 #undef __FUNCT__
885 #define __FUNCT__ "SNESLineSearchGetSNES"
886 /*@
887    SNESLineSearchGetSNES - Gets the SNES instance associated with the line search.
888    Having an associated SNES is necessary because most line search implementations must be able to
889    evaluate the function using SNESComputeFunction() for the associated SNES.  This routine
890    is used in the line search implementations when one must get this associated SNES instance.
891 
892    Input Parameters:
893 .  linesearch - linesearch context
894 
895    Output Parameters:
896 .  snes - The snes instance
897 
898    Level: developer
899 
900 .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
901 @*/
902 PetscErrorCode  SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes)
903 {
904   PetscFunctionBegin;
905   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
906   PetscValidPointer(snes, 2);
907   *snes = linesearch->snes;
908   PetscFunctionReturn(0);
909 }
910 
911 #undef __FUNCT__
912 #define __FUNCT__ "SNESLineSearchGetLambda"
913 /*@
914    SNESLineSearchGetLambda - Gets the last linesearch steplength discovered.
915 
916    Input Parameters:
917 .  linesearch - linesearch context
918 
919    Output Parameters:
920 .  lambda - The last steplength computed during SNESLineSearchApply()
921 
922    Level: advanced
923 
924    Notes:
925    This is useful in methods where the solver is ill-scaled and
926    requires some adaptive notion of the difference in scale between the
927    solution and the function.  For instance, SNESQN may be scaled by the
928    line search lambda using the argument -snes_qn_scaling ls.
929 
930 
931 .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetDamping(), SNESLineSearchApply()
932 @*/
933 PetscErrorCode  SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda)
934 {
935   PetscFunctionBegin;
936   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
937   PetscValidPointer(lambda, 2);
938   *lambda = linesearch->lambda;
939   PetscFunctionReturn(0);
940 }
941 
942 #undef __FUNCT__
943 #define __FUNCT__ "SNESLineSearchSetLambda"
944 /*@
945    SNESLineSearchSetLambda - Sets the linesearch steplength.
946 
947    Input Parameters:
948 +  linesearch - linesearch context
949 -  lambda - The last steplength.
950 
951    Notes:
952    This routine is typically used within implementations of SNESLineSearchApply
953    to set the final steplength.  This routine (and SNESLineSearchGetLambda()) were
954    added in order to facilitate Quasi-Newton methods that use the previous steplength
955    as an inner scaling parameter.
956 
957    Level: advanced
958 
959 .seealso: SNESLineSearchGetLambda()
960 @*/
961 PetscErrorCode  SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
962 {
963   PetscFunctionBegin;
964   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
965   linesearch->lambda = lambda;
966   PetscFunctionReturn(0);
967 }
968 
969 #undef  __FUNCT__
970 #define __FUNCT__ "SNESLineSearchGetTolerances"
971 /*@
972    SNESLineSearchGetTolerances - Gets the tolerances for the linesearch.  These include
973    tolerances for the relative and absolute change in the function norm, the change
974    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
975    and the maximum number of iterations the line search procedure may take.
976 
977    Input Parameters:
978 .  linesearch - linesearch context
979 
980    Output Parameters:
981 +  steptol - The minimum steplength
982 .  maxstep - The maximum steplength
983 .  rtol    - The relative tolerance for iterative line searches
984 .  atol    - The absolute tolerance for iterative line searches
985 .  ltol    - The change in lambda tolerance for iterative line searches
986 -  max_it  - The maximum number of iterations of the line search
987 
988    Level: intermediate
989 
990    Notes:
991    Different line searches may implement these parameters slightly differently as
992    the type requires.
993 
994 .seealso: SNESLineSearchSetTolerances()
995 @*/
996 PetscErrorCode  SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
997 {
998   PetscFunctionBegin;
999   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1000   if (steptol) {
1001     PetscValidPointer(steptol, 2);
1002     *steptol = linesearch->steptol;
1003   }
1004   if (maxstep) {
1005     PetscValidPointer(maxstep, 3);
1006     *maxstep = linesearch->maxstep;
1007   }
1008   if (rtol) {
1009     PetscValidPointer(rtol, 4);
1010     *rtol = linesearch->rtol;
1011   }
1012   if (atol) {
1013     PetscValidPointer(atol, 5);
1014     *atol = linesearch->atol;
1015   }
1016   if (ltol) {
1017     PetscValidPointer(ltol, 6);
1018     *ltol = linesearch->ltol;
1019   }
1020   if (max_its) {
1021     PetscValidPointer(max_its, 7);
1022     *max_its = linesearch->max_its;
1023   }
1024   PetscFunctionReturn(0);
1025 }
1026 
1027 #undef  __FUNCT__
1028 #define __FUNCT__ "SNESLineSearchSetTolerances"
1029 /*@
1030    SNESLineSearchSetTolerances -  Gets the tolerances for the linesearch.  These include
1031    tolerances for the relative and absolute change in the function norm, the change
1032    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
1033    and the maximum number of iterations the line search procedure may take.
1034 
1035    Input Parameters:
1036 +  linesearch - linesearch context
1037 .  steptol - The minimum steplength
1038 .  maxstep - The maximum steplength
1039 .  rtol    - The relative tolerance for iterative line searches
1040 .  atol    - The absolute tolerance for iterative line searches
1041 .  ltol    - The change in lambda tolerance for iterative line searches
1042 -  max_it  - The maximum number of iterations of the line search
1043 
1044    Notes:
1045    The user may choose to not set any of the tolerances using PETSC_DEFAULT in
1046    place of an argument.
1047 
1048    Level: intermediate
1049 
1050 .seealso: SNESLineSearchGetTolerances()
1051 @*/
1052 PetscErrorCode  SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its)
1053 {
1054   PetscFunctionBegin;
1055   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1056   PetscValidLogicalCollectiveReal(linesearch,steptol,2);
1057   PetscValidLogicalCollectiveReal(linesearch,maxstep,3);
1058   PetscValidLogicalCollectiveReal(linesearch,rtol,4);
1059   PetscValidLogicalCollectiveReal(linesearch,atol,5);
1060   PetscValidLogicalCollectiveReal(linesearch,ltol,6);
1061   PetscValidLogicalCollectiveInt(linesearch,max_its,7);
1062 
1063   if (steptol!= PETSC_DEFAULT) {
1064     if (steptol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Minimum step length %14.12e must be non-negative",(double)steptol);
1065     linesearch->steptol = steptol;
1066   }
1067 
1068   if (maxstep!= PETSC_DEFAULT) {
1069     if (maxstep < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum step length %14.12e must be non-negative",(double)maxstep);
1070     linesearch->maxstep = maxstep;
1071   }
1072 
1073   if (rtol != PETSC_DEFAULT) {
1074     if (rtol < 0.0 || 1.0 <= rtol) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Relative tolerance %14.12e must be non-negative and less than 1.0",(double)rtol);
1075     linesearch->rtol = rtol;
1076   }
1077 
1078   if (atol != PETSC_DEFAULT) {
1079     if (atol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %14.12e must be non-negative",(double)atol);
1080     linesearch->atol = atol;
1081   }
1082 
1083   if (ltol != PETSC_DEFAULT) {
1084     if (ltol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Labmda tolerance %14.12e must be non-negative",(double)ltol);
1085     linesearch->ltol = ltol;
1086   }
1087 
1088   if (max_its != PETSC_DEFAULT) {
1089     if (max_its < 0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",max_its);
1090     linesearch->max_its = max_its;
1091   }
1092   PetscFunctionReturn(0);
1093 }
1094 
1095 #undef __FUNCT__
1096 #define __FUNCT__ "SNESLineSearchGetDamping"
1097 /*@
1098    SNESLineSearchGetDamping - Gets the line search damping parameter.
1099 
1100    Input Parameters:
1101 .  linesearch - linesearch context
1102 
1103    Output Parameters:
1104 .  damping - The damping parameter
1105 
1106    Level: advanced
1107 
1108 .seealso: SNESLineSearchGetStepTolerance(), SNESQN
1109 @*/
1110 
1111 PetscErrorCode  SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping)
1112 {
1113   PetscFunctionBegin;
1114   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1115   PetscValidPointer(damping, 2);
1116   *damping = linesearch->damping;
1117   PetscFunctionReturn(0);
1118 }
1119 
1120 #undef __FUNCT__
1121 #define __FUNCT__ "SNESLineSearchSetDamping"
1122 /*@
1123    SNESLineSearchSetDamping - Sets the line search damping paramter.
1124 
1125    Input Parameters:
1126 .  linesearch - linesearch context
1127 .  damping - The damping parameter
1128 
1129    Level: intermediate
1130 
1131    Notes:
1132    The basic line search merely takes the update step scaled by the damping parameter.
1133    The use of the damping parameter in the l2 and cp line searches is much more subtle;
1134    it is used as a starting point in calculating the secant step. However, the eventual
1135    step may be of greater length than the damping parameter.  In the bt line search it is
1136    used as the maximum possible step length, as the bt line search only backtracks.
1137 
1138 .seealso: SNESLineSearchGetDamping()
1139 @*/
1140 PetscErrorCode  SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping)
1141 {
1142   PetscFunctionBegin;
1143   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1144   linesearch->damping = damping;
1145   PetscFunctionReturn(0);
1146 }
1147 
1148 #undef __FUNCT__
1149 #define __FUNCT__ "SNESLineSearchGetOrder"
1150 /*@
1151    SNESLineSearchGetOrder - Gets the line search approximation order.
1152 
1153    Input Parameters:
1154 .  linesearch - linesearch context
1155 
1156    Output Parameters:
1157 .  order - The order
1158 
1159    Possible Values for order:
1160 +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
1161 .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
1162 -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
1163 
1164    Level: intermediate
1165 
1166 .seealso: SNESLineSearchSetOrder()
1167 @*/
1168 
1169 PetscErrorCode  SNESLineSearchGetOrder(SNESLineSearch linesearch,PetscInt *order)
1170 {
1171   PetscFunctionBegin;
1172   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1173   PetscValidPointer(order, 2);
1174   *order = linesearch->order;
1175   PetscFunctionReturn(0);
1176 }
1177 
1178 #undef __FUNCT__
1179 #define __FUNCT__ "SNESLineSearchSetOrder"
1180 /*@
1181    SNESLineSearchSetOrder - Sets the line search damping paramter.
1182 
1183    Input Parameters:
1184 .  linesearch - linesearch context
1185 .  order - The damping parameter
1186 
1187    Level: intermediate
1188 
1189    Possible Values for order:
1190 +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
1191 .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
1192 -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
1193 
1194    Notes:
1195    Variable orders are supported by the following line searches:
1196 +  bt - cubic and quadratic
1197 -  cp - linear and quadratic
1198 
1199 .seealso: SNESLineSearchGetOrder()
1200 @*/
1201 PetscErrorCode  SNESLineSearchSetOrder(SNESLineSearch linesearch,PetscInt order)
1202 {
1203   PetscFunctionBegin;
1204   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1205   linesearch->order = order;
1206   PetscFunctionReturn(0);
1207 }
1208 
1209 #undef __FUNCT__
1210 #define __FUNCT__ "SNESLineSearchGetNorms"
1211 /*@
1212    SNESLineSearchGetNorms - Gets the norms for for X, Y, and F.
1213 
1214    Input Parameters:
1215 .  linesearch - linesearch context
1216 
1217    Output Parameters:
1218 +  xnorm - The norm of the current solution
1219 .  fnorm - The norm of the current function
1220 -  ynorm - The norm of the current update
1221 
1222    Notes:
1223    This function is mainly called from SNES implementations.
1224 
1225    Level: developer
1226 
1227 .seealso: SNESLineSearchSetNorms() SNESLineSearchGetVecs()
1228 @*/
1229 PetscErrorCode  SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm)
1230 {
1231   PetscFunctionBegin;
1232   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1233   if (xnorm) *xnorm = linesearch->xnorm;
1234   if (fnorm) *fnorm = linesearch->fnorm;
1235   if (ynorm) *ynorm = linesearch->ynorm;
1236   PetscFunctionReturn(0);
1237 }
1238 
1239 #undef __FUNCT__
1240 #define __FUNCT__ "SNESLineSearchSetNorms"
1241 /*@
1242    SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F.
1243 
1244    Input Parameters:
1245 +  linesearch - linesearch context
1246 .  xnorm - The norm of the current solution
1247 .  fnorm - The norm of the current function
1248 -  ynorm - The norm of the current update
1249 
1250    Level: advanced
1251 
1252 .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1253 @*/
1254 PetscErrorCode  SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
1255 {
1256   PetscFunctionBegin;
1257   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1258   linesearch->xnorm = xnorm;
1259   linesearch->fnorm = fnorm;
1260   linesearch->ynorm = ynorm;
1261   PetscFunctionReturn(0);
1262 }
1263 
1264 #undef __FUNCT__
1265 #define __FUNCT__ "SNESLineSearchComputeNorms"
1266 /*@
1267    SNESLineSearchComputeNorms - Computes the norms of X, F, and Y.
1268 
1269    Input Parameters:
1270 .  linesearch - linesearch context
1271 
1272    Options Database Keys:
1273 .   -snes_linesearch_norms - turn norm computation on or off
1274 
1275    Level: intermediate
1276 
1277 .seealso: SNESLineSearchGetNorms, SNESLineSearchSetNorms(), SNESLineSearchSetComputeNorms()
1278 @*/
1279 PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
1280 {
1281   PetscErrorCode ierr;
1282   SNES           snes;
1283 
1284   PetscFunctionBegin;
1285   if (linesearch->norms) {
1286     if (linesearch->ops->vinorm) {
1287       ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
1288       ierr = VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
1289       ierr = VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
1290       ierr = (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);CHKERRQ(ierr);
1291     } else {
1292       ierr = VecNormBegin(linesearch->vec_func,   NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
1293       ierr = VecNormBegin(linesearch->vec_sol,    NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
1294       ierr = VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
1295       ierr = VecNormEnd(linesearch->vec_func,     NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
1296       ierr = VecNormEnd(linesearch->vec_sol,      NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
1297       ierr = VecNormEnd(linesearch->vec_update,   NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
1298     }
1299   }
1300   PetscFunctionReturn(0);
1301 }
1302 
1303 #undef __FUNCT__
1304 #define __FUNCT__ "SNESLineSearchSetComputeNorms"
1305 /*@
1306    SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search.
1307 
1308    Input Parameters:
1309 +  linesearch  - linesearch context
1310 -  flg  - indicates whether or not to compute norms
1311 
1312    Options Database Keys:
1313 .   -snes_linesearch_norms - turn norm computation on or off
1314 
1315    Notes:
1316    This is most relevant to the SNESLINESEARCHBASIC line search type.
1317 
1318    Level: intermediate
1319 
1320 .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetNorms(), SNESLineSearchComputeNorms(), SNESLINESEARCHBASIC
1321 @*/
1322 PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg)
1323 {
1324   PetscFunctionBegin;
1325   linesearch->norms = flg;
1326   PetscFunctionReturn(0);
1327 }
1328 
1329 #undef __FUNCT__
1330 #define __FUNCT__ "SNESLineSearchGetVecs"
1331 /*@
1332    SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context
1333 
1334    Input Parameters:
1335 .  linesearch - linesearch context
1336 
1337    Output Parameters:
1338 +  X - Solution vector
1339 .  F - Function vector
1340 .  Y - Search direction vector
1341 .  W - Solution work vector
1342 -  G - Function work vector
1343 
1344    Notes:
1345    At the beginning of a line search application, X should contain a
1346    solution and the vector F the function computed at X.  At the end of the
1347    line search application, X should contain the new solution, and F the
1348    function evaluated at the new solution.
1349 
1350    Level: advanced
1351 
1352 .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1353 @*/
1354 PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G)
1355 {
1356   PetscFunctionBegin;
1357   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1358   if (X) {
1359     PetscValidPointer(X, 2);
1360     *X = linesearch->vec_sol;
1361   }
1362   if (F) {
1363     PetscValidPointer(F, 3);
1364     *F = linesearch->vec_func;
1365   }
1366   if (Y) {
1367     PetscValidPointer(Y, 4);
1368     *Y = linesearch->vec_update;
1369   }
1370   if (W) {
1371     PetscValidPointer(W, 5);
1372     *W = linesearch->vec_sol_new;
1373   }
1374   if (G) {
1375     PetscValidPointer(G, 6);
1376     *G = linesearch->vec_func_new;
1377   }
1378   PetscFunctionReturn(0);
1379 }
1380 
1381 #undef __FUNCT__
1382 #define __FUNCT__ "SNESLineSearchSetVecs"
1383 /*@
1384    SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context
1385 
1386    Input Parameters:
1387 +  linesearch - linesearch context
1388 .  X - Solution vector
1389 .  F - Function vector
1390 .  Y - Search direction vector
1391 .  W - Solution work vector
1392 -  G - Function work vector
1393 
1394    Level: advanced
1395 
1396 .seealso: SNESLineSearchSetNorms(), SNESLineSearchGetVecs()
1397 @*/
1398 PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G)
1399 {
1400   PetscFunctionBegin;
1401   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1402   if (X) {
1403     PetscValidHeaderSpecific(X,VEC_CLASSID,2);
1404     linesearch->vec_sol = X;
1405   }
1406   if (F) {
1407     PetscValidHeaderSpecific(F,VEC_CLASSID,3);
1408     linesearch->vec_func = F;
1409   }
1410   if (Y) {
1411     PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
1412     linesearch->vec_update = Y;
1413   }
1414   if (W) {
1415     PetscValidHeaderSpecific(W,VEC_CLASSID,5);
1416     linesearch->vec_sol_new = W;
1417   }
1418   if (G) {
1419     PetscValidHeaderSpecific(G,VEC_CLASSID,6);
1420     linesearch->vec_func_new = G;
1421   }
1422   PetscFunctionReturn(0);
1423 }
1424 
1425 #undef __FUNCT__
1426 #define __FUNCT__ "SNESLineSearchAppendOptionsPrefix"
1427 /*@C
1428    SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1429    SNES options in the database.
1430 
1431    Logically Collective on SNESLineSearch
1432 
1433    Input Parameters:
1434 +  snes - the SNES context
1435 -  prefix - the prefix to prepend to all option names
1436 
1437    Notes:
1438    A hyphen (-) must NOT be given at the beginning of the prefix name.
1439    The first character of all runtime options is AUTOMATICALLY the hyphen.
1440 
1441    Level: advanced
1442 
1443 .keywords: SNESLineSearch, append, options, prefix, database
1444 
1445 .seealso: SNESGetOptionsPrefix()
1446 @*/
1447 PetscErrorCode  SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[])
1448 {
1449   PetscErrorCode ierr;
1450 
1451   PetscFunctionBegin;
1452   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1453   ierr = PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1454   PetscFunctionReturn(0);
1455 }
1456 
1457 #undef __FUNCT__
1458 #define __FUNCT__ "SNESLineSearchGetOptionsPrefix"
1459 /*@C
1460    SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all
1461    SNESLineSearch options in the database.
1462 
1463    Not Collective
1464 
1465    Input Parameter:
1466 .  linesearch - the SNESLineSearch context
1467 
1468    Output Parameter:
1469 .  prefix - pointer to the prefix string used
1470 
1471    Notes:
1472    On the fortran side, the user should pass in a string 'prefix' of
1473    sufficient length to hold the prefix.
1474 
1475    Level: advanced
1476 
1477 .keywords: SNESLineSearch, get, options, prefix, database
1478 
1479 .seealso: SNESAppendOptionsPrefix()
1480 @*/
1481 PetscErrorCode  SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[])
1482 {
1483   PetscErrorCode ierr;
1484 
1485   PetscFunctionBegin;
1486   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1487   ierr = PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1488   PetscFunctionReturn(0);
1489 }
1490 
1491 #undef __FUNCT__
1492 #define __FUNCT__ "SNESLineSearchSetWorkVecs"
1493 /*@C
1494    SNESLineSearchSetWorkVecs - Gets work vectors for the line search.
1495 
1496    Input Parameter:
1497 +  linesearch - the SNESLineSearch context
1498 -  nwork - the number of work vectors
1499 
1500    Level: developer
1501 
1502    Developers Note: This is PETSC_EXTERN because it may be used by user written plugin SNES implementations
1503 
1504 .keywords: SNESLineSearch, work, vector
1505 
1506 .seealso: SNESSetWorkVecs()
1507 @*/
1508 PetscErrorCode  SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork)
1509 {
1510   PetscErrorCode ierr;
1511 
1512   PetscFunctionBegin;
1513   if (linesearch->vec_sol) {
1514     ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr);
1515   } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
1516   PetscFunctionReturn(0);
1517 }
1518 
1519 #undef __FUNCT__
1520 #define __FUNCT__ "SNESLineSearchGetSuccess"
1521 /*@
1522    SNESLineSearchGetSuccess - Gets the success/failure status of the last line search application
1523 
1524    Input Parameters:
1525 .  linesearch - linesearch context
1526 
1527    Output Parameters:
1528 .  success - The success or failure status
1529 
1530    Notes:
1531    This is typically called after SNESLineSearchApply() in order to determine if the line-search failed
1532    (and set the SNES convergence accordingly).
1533 
1534    Level: intermediate
1535 
1536 .seealso: SNESLineSearchSetSuccess()
1537 @*/
1538 PetscErrorCode  SNESLineSearchGetSuccess(SNESLineSearch linesearch, PetscBool *success)
1539 {
1540   PetscFunctionBegin;
1541   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1542   PetscValidPointer(success, 2);
1543   if (success) *success = linesearch->success;
1544   PetscFunctionReturn(0);
1545 }
1546 
1547 #undef __FUNCT__
1548 #define __FUNCT__ "SNESLineSearchSetSuccess"
1549 /*@
1550    SNESLineSearchSetSuccess - Sets the success/failure status of the last line search application
1551 
1552    Input Parameters:
1553 +  linesearch - linesearch context
1554 -  success - The success or failure status
1555 
1556    Notes:
1557    This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set
1558    the success or failure of the line search method.
1559 
1560    Level: developer
1561 
1562 .seealso: SNESLineSearchGetSuccess()
1563 @*/
1564 PetscErrorCode  SNESLineSearchSetSuccess(SNESLineSearch linesearch, PetscBool success)
1565 {
1566   PetscFunctionBegin;
1567   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1568   linesearch->success = success;
1569   PetscFunctionReturn(0);
1570 }
1571 
1572 #undef __FUNCT__
1573 #define __FUNCT__ "SNESLineSearchSetVIFunctions"
1574 /*@C
1575    SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.
1576 
1577    Input Parameters:
1578 +  snes - nonlinear context obtained from SNESCreate()
1579 .  projectfunc - function for projecting the function to the bounds
1580 -  normfunc - function for computing the norm of an active set
1581 
1582    Logically Collective on SNES
1583 
1584    Calling sequence of projectfunc:
1585 .vb
1586    projectfunc (SNES snes, Vec X)
1587 .ve
1588 
1589     Input parameters for projectfunc:
1590 +   snes - nonlinear context
1591 -   X - current solution
1592 
1593     Output parameters for projectfunc:
1594 .   X - Projected solution
1595 
1596    Calling sequence of normfunc:
1597 .vb
1598    projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm)
1599 .ve
1600 
1601     Input parameters for normfunc:
1602 +   snes - nonlinear context
1603 .   X - current solution
1604 -   F - current residual
1605 
1606     Output parameters for normfunc:
1607 .   fnorm - VI-specific norm of the function
1608 
1609     Notes:
1610     The VI solvers require projection of the solution to the feasible set.  projectfunc should implement this.
1611 
1612     The VI solvers require special evaluation of the function norm such that the norm is only calculated
1613     on the inactive set.  This should be implemented by normfunc.
1614 
1615     Level: developer
1616 
1617 .keywords: SNES, line search, VI, nonlinear, set, line search
1618 
1619 .seealso: SNESLineSearchGetVIFunctions(), SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck()
1620 @*/
1621 extern PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc)
1622 {
1623   PetscFunctionBegin;
1624   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1625   if (projectfunc) linesearch->ops->viproject = projectfunc;
1626   if (normfunc) linesearch->ops->vinorm = normfunc;
1627   PetscFunctionReturn(0);
1628 }
1629 
1630 #undef __FUNCT__
1631 #define __FUNCT__ "SNESLineSearchGetVIFunctions"
1632 /*@C
1633    SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.
1634 
1635    Input Parameters:
1636 .  snes - nonlinear context obtained from SNESCreate()
1637 
1638    Output Parameters:
1639 +  projectfunc - function for projecting the function to the bounds
1640 -  normfunc - function for computing the norm of an active set
1641 
1642    Logically Collective on SNES
1643 
1644     Level: developer
1645 
1646 .keywords: SNES, line search, VI, nonlinear, get, line search
1647 
1648 .seealso: SNESLineSearchSetVIFunctions(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
1649 @*/
1650 extern PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc)
1651 {
1652   PetscFunctionBegin;
1653   if (projectfunc) *projectfunc = linesearch->ops->viproject;
1654   if (normfunc) *normfunc = linesearch->ops->vinorm;
1655   PetscFunctionReturn(0);
1656 }
1657 
1658 #undef __FUNCT__
1659 #define __FUNCT__ "SNESLineSearchRegister"
1660 /*@C
1661   SNESLineSearchRegister - See SNESLineSearchRegister()
1662 
1663   Level: advanced
1664 @*/
1665 PetscErrorCode  SNESLineSearchRegister(const char sname[],PetscErrorCode (*function)(SNESLineSearch))
1666 {
1667   PetscErrorCode ierr;
1668 
1669   PetscFunctionBegin;
1670   ierr = PetscFunctionListAdd(&SNESLineSearchList,sname,function);CHKERRQ(ierr);
1671   PetscFunctionReturn(0);
1672 }
1673