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