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