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