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