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