xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision c197fec04cad292fef96e4490994e90a01637e61)
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->maxstep      = 1e8;
196   linesearch->steptol      = 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_its      = 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, l2, 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 search length 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 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`, `l2`, `cp`, `nleqerr`, `bisection`, `shell`
797 . -snes_linesearch_order <order>                                    - 1, 2, 3.  Most types only support certain orders (`bt` supports 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 step length
800 . -snes_linesearch_maxstep                                          - The maximum step size
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 search length 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 step length", "SNESLineSearchSetTolerances", linesearch->steptol, &linesearch->steptol, NULL));
842   PetscCall(PetscOptionsReal("-snes_linesearch_maxstep", "Maximum step size", "SNESLineSearchSetTolerances", linesearch->maxstep, &linesearch->maxstep, 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_its, &linesearch->max_its, NULL));
847 
848   /* damping parameters */
849   PetscCall(PetscOptionsReal("-snes_linesearch_damping", "Line search damping and initial step guess", "SNESLineSearchSetDamping", linesearch->damping, &linesearch->damping, NULL));
850 
851   PetscCall(PetscOptionsBool("-snes_linesearch_keeplambda", "Use previous lambda as damping", "SNESLineSearchSetKeepLambda", linesearch->keeplambda, &linesearch->keeplambda, NULL));
852 
853   /* precheck */
854   PetscCall(PetscOptionsBool("-snes_linesearch_precheck_picard", "Use a correction that sometimes improves convergence of Picard iteration", "SNESLineSearchPreCheckPicard", flg, &flg, &set));
855   if (set) {
856     if (flg) {
857       linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */
858 
859       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));
860       PetscCall(SNESLineSearchSetPreCheck(linesearch, SNESLineSearchPreCheckPicard, &linesearch->precheck_picard_angle));
861     } else {
862       PetscCall(SNESLineSearchSetPreCheck(linesearch, NULL, NULL));
863     }
864   }
865   PetscCall(PetscOptionsInt("-snes_linesearch_order", "Order of approximation used in the line search", "SNESLineSearchSetOrder", linesearch->order, &linesearch->order, NULL));
866   PetscCall(PetscOptionsBool("-snes_linesearch_norms", "Compute final norms in line search", "SNESLineSearchSetComputeNorms", linesearch->norms, &linesearch->norms, NULL));
867 
868   PetscTryTypeMethod(linesearch, setfromoptions, PetscOptionsObject);
869 
870   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)linesearch, PetscOptionsObject));
871   PetscOptionsEnd();
872   PetscFunctionReturn(PETSC_SUCCESS);
873 }
874 
875 /*@
876   SNESLineSearchView - Prints useful information about the line search
877 
878   Logically Collective
879 
880   Input Parameters:
881 + linesearch - line search context
882 - viewer     - the `PetscViewer` to display the line search information to
883 
884   Level: intermediate
885 
886 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `PetscViewer`, `SNESLineSearchCreate()`
887 @*/
888 PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer)
889 {
890   PetscBool iascii;
891 
892   PetscFunctionBegin;
893   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
894   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch), &viewer));
895   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
896   PetscCheckSameComm(linesearch, 1, viewer, 2);
897 
898   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
899   if (iascii) {
900     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)linesearch, viewer));
901     PetscCall(PetscViewerASCIIPushTab(viewer));
902     PetscTryTypeMethod(linesearch, view, viewer);
903     PetscCall(PetscViewerASCIIPopTab(viewer));
904     PetscCall(PetscViewerASCIIPrintf(viewer, "  maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep, (double)linesearch->steptol));
905     PetscCall(PetscViewerASCIIPrintf(viewer, "  tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol, (double)linesearch->atol, (double)linesearch->ltol));
906     PetscCall(PetscViewerASCIIPrintf(viewer, "  maximum iterations=%" PetscInt_FMT "\n", linesearch->max_its));
907     if (linesearch->ops->precheck) {
908       if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) {
909         PetscCall(PetscViewerASCIIPrintf(viewer, "  using precheck step to speed up Picard convergence\n"));
910       } else {
911         PetscCall(PetscViewerASCIIPrintf(viewer, "  using user-defined precheck step\n"));
912       }
913     }
914     if (linesearch->ops->postcheck) PetscCall(PetscViewerASCIIPrintf(viewer, "  using user-defined postcheck step\n"));
915   }
916   PetscFunctionReturn(PETSC_SUCCESS);
917 }
918 
919 /*@
920   SNESLineSearchGetType - Gets the `SNESLinesearchType` of a `SNESLineSearch`
921 
922   Logically Collective
923 
924   Input Parameter:
925 . linesearch - the line search context
926 
927   Output Parameter:
928 . type - The type of line search, or `NULL` if not set
929 
930   Level: intermediate
931 
932 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchSetType()`
933 @*/
934 PetscErrorCode SNESLineSearchGetType(SNESLineSearch linesearch, SNESLineSearchType *type)
935 {
936   PetscFunctionBegin;
937   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
938   PetscAssertPointer(type, 2);
939   *type = ((PetscObject)linesearch)->type_name;
940   PetscFunctionReturn(PETSC_SUCCESS);
941 }
942 
943 /*@
944   SNESLineSearchSetType - Sets the `SNESLinesearchType` of a `SNESLineSearch` object to indicate the line search algorithm that should be used by a given `SNES` solver
945 
946   Logically Collective
947 
948   Input Parameters:
949 + linesearch - the line search context
950 - type       - The type of line search to be used, see `SNESLineSearchType`
951 
952   Options Database Key:
953 . -snes_linesearch_type <type> - basic (or equivalently none), bt, l2, cp, nleqerr, shell
954 
955   Level: intermediate
956 
957   Note:
958   The `SNESLineSearch` object is generally obtained with `SNESGetLineSearch()`
959 
960 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchGetType()`,
961           `SNESGetLineSearch()`
962 @*/
963 PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type)
964 {
965   PetscBool match;
966   PetscErrorCode (*r)(SNESLineSearch);
967 
968   PetscFunctionBegin;
969   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
970   PetscAssertPointer(type, 2);
971 
972   PetscCall(PetscObjectTypeCompare((PetscObject)linesearch, type, &match));
973   if (match) PetscFunctionReturn(PETSC_SUCCESS);
974 
975   PetscCall(PetscFunctionListFind(SNESLineSearchList, type, &r));
976   PetscCheck(r, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested Line Search type %s", type);
977   /* Destroy the previous private line search context */
978   PetscTryTypeMethod(linesearch, destroy);
979   linesearch->ops->destroy = NULL;
980   /* Reinitialize function pointers in SNESLineSearchOps structure */
981   linesearch->ops->apply          = NULL;
982   linesearch->ops->view           = NULL;
983   linesearch->ops->setfromoptions = NULL;
984   linesearch->ops->destroy        = NULL;
985 
986   PetscCall(PetscObjectChangeTypeName((PetscObject)linesearch, type));
987   PetscCall((*r)(linesearch));
988   PetscFunctionReturn(PETSC_SUCCESS);
989 }
990 
991 /*@
992   SNESLineSearchSetSNES - Sets the `SNES` for the linesearch for function evaluation.
993 
994   Input Parameters:
995 + linesearch - the line search context
996 - snes       - The `SNES` instance
997 
998   Level: developer
999 
1000   Note:
1001   This happens automatically when the line search is obtained/created with
1002   `SNESGetLineSearch()`.  This routine is therefore mainly called within `SNES`
1003   implementations.
1004 
1005 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetSNES()`, `SNESLineSearchSetVecs()`
1006 @*/
1007 PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes)
1008 {
1009   PetscFunctionBegin;
1010   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1011   PetscValidHeaderSpecific(snes, SNES_CLASSID, 2);
1012   linesearch->snes = snes;
1013   PetscFunctionReturn(PETSC_SUCCESS);
1014 }
1015 
1016 /*@
1017   SNESLineSearchGetSNES - Gets the `SNES` instance associated with the line search.
1018 
1019   Not Collective
1020 
1021   Input Parameter:
1022 . linesearch - the line search context
1023 
1024   Output Parameter:
1025 . snes - The `SNES` instance
1026 
1027   Level: developer
1028 
1029 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESType`, `SNESLineSearchSetVecs()`
1030 @*/
1031 PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes)
1032 {
1033   PetscFunctionBegin;
1034   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1035   PetscAssertPointer(snes, 2);
1036   *snes = linesearch->snes;
1037   PetscFunctionReturn(PETSC_SUCCESS);
1038 }
1039 
1040 /*@
1041   SNESLineSearchGetLambda - Gets the last line search steplength used
1042 
1043   Not Collective
1044 
1045   Input Parameter:
1046 . linesearch - the line search context
1047 
1048   Output Parameter:
1049 . lambda - The last steplength computed during `SNESLineSearchApply()`
1050 
1051   Level: advanced
1052 
1053   Note:
1054   This is useful in methods where the solver is ill-scaled and
1055   requires some adaptive notion of the difference in scale between the
1056   solution and the function.  For instance, `SNESQN` may be scaled by the
1057   line search lambda using the argument -snes_qn_scaling ls.
1058 
1059 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetLambda()`, `SNESLineSearchGetDamping()`, `SNESLineSearchApply()`
1060 @*/
1061 PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch linesearch, PetscReal *lambda)
1062 {
1063   PetscFunctionBegin;
1064   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1065   PetscAssertPointer(lambda, 2);
1066   *lambda = linesearch->lambda;
1067   PetscFunctionReturn(PETSC_SUCCESS);
1068 }
1069 
1070 /*@
1071   SNESLineSearchSetLambda - Sets the line search steplength
1072 
1073   Input Parameters:
1074 + linesearch - line search context
1075 - lambda     - The steplength to use
1076 
1077   Level: advanced
1078 
1079   Note:
1080   This routine is typically used within implementations of `SNESLineSearchApply()`
1081   to set the final steplength.  This routine (and `SNESLineSearchGetLambda()`) were
1082   added in order to facilitate Quasi-Newton methods that use the previous steplength
1083   as an inner scaling parameter.
1084 
1085 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetLambda()`
1086 @*/
1087 PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
1088 {
1089   PetscFunctionBegin;
1090   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1091   linesearch->lambda = lambda;
1092   PetscFunctionReturn(PETSC_SUCCESS);
1093 }
1094 
1095 /*@
1096   SNESLineSearchGetTolerances - Gets the tolerances for the line search.
1097 
1098   Not Collective
1099 
1100   Input Parameter:
1101 . linesearch - the line search context
1102 
1103   Output Parameters:
1104 + steptol - The minimum steplength
1105 . maxstep - The maximum steplength
1106 . rtol    - The relative tolerance for iterative line searches
1107 . atol    - The absolute tolerance for iterative line searches
1108 . ltol    - The change in lambda tolerance for iterative line searches
1109 - max_its - The maximum number of iterations of the line search
1110 
1111   Level: intermediate
1112 
1113   Note:
1114   Different line searches may implement these parameters slightly differently as
1115   the type requires.
1116 
1117 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetTolerances()`
1118 @*/
1119 PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch linesearch, PetscReal *steptol, PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
1120 {
1121   PetscFunctionBegin;
1122   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1123   if (steptol) {
1124     PetscAssertPointer(steptol, 2);
1125     *steptol = linesearch->steptol;
1126   }
1127   if (maxstep) {
1128     PetscAssertPointer(maxstep, 3);
1129     *maxstep = linesearch->maxstep;
1130   }
1131   if (rtol) {
1132     PetscAssertPointer(rtol, 4);
1133     *rtol = linesearch->rtol;
1134   }
1135   if (atol) {
1136     PetscAssertPointer(atol, 5);
1137     *atol = linesearch->atol;
1138   }
1139   if (ltol) {
1140     PetscAssertPointer(ltol, 6);
1141     *ltol = linesearch->ltol;
1142   }
1143   if (max_its) {
1144     PetscAssertPointer(max_its, 7);
1145     *max_its = linesearch->max_its;
1146   }
1147   PetscFunctionReturn(PETSC_SUCCESS);
1148 }
1149 
1150 /*@
1151   SNESLineSearchSetTolerances -  Sets the tolerances for the linesearch.
1152 
1153   Collective
1154 
1155   Input Parameters:
1156 + linesearch - the line search context
1157 . steptol    - The minimum steplength
1158 . maxstep    - The maximum steplength
1159 . rtol       - The relative tolerance for iterative line searches
1160 . atol       - The absolute tolerance for iterative line searches
1161 . ltol       - The change in lambda tolerance for iterative line searches
1162 - max_it     - The maximum number of iterations of the line search
1163 
1164   Options Database Keys:
1165 + -snes_linesearch_minlambda - The minimum step length
1166 . -snes_linesearch_maxstep   - The maximum step size
1167 . -snes_linesearch_rtol      - Relative tolerance for iterative line searches
1168 . -snes_linesearch_atol      - Absolute tolerance for iterative line searches
1169 . -snes_linesearch_ltol      - Change in lambda tolerance for iterative line searches
1170 - -snes_linesearch_max_it    - The number of iterations for iterative line searches
1171 
1172   Level: intermediate
1173 
1174   Note:
1175   The user may choose to not set any of the tolerances using `PETSC_DEFAULT` in place of an argument.
1176 
1177 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetTolerances()`
1178 @*/
1179 PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch linesearch, PetscReal steptol, PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_it)
1180 {
1181   PetscFunctionBegin;
1182   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1183   PetscValidLogicalCollectiveReal(linesearch, steptol, 2);
1184   PetscValidLogicalCollectiveReal(linesearch, maxstep, 3);
1185   PetscValidLogicalCollectiveReal(linesearch, rtol, 4);
1186   PetscValidLogicalCollectiveReal(linesearch, atol, 5);
1187   PetscValidLogicalCollectiveReal(linesearch, ltol, 6);
1188   PetscValidLogicalCollectiveInt(linesearch, max_it, 7);
1189 
1190   if (steptol != (PetscReal)PETSC_DEFAULT) {
1191     PetscCheck(steptol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Minimum step length %14.12e must be non-negative", (double)steptol);
1192     linesearch->steptol = steptol;
1193   }
1194 
1195   if (maxstep != (PetscReal)PETSC_DEFAULT) {
1196     PetscCheck(maxstep >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum step length %14.12e must be non-negative", (double)maxstep);
1197     linesearch->maxstep = maxstep;
1198   }
1199 
1200   if (rtol != (PetscReal)PETSC_DEFAULT) {
1201     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);
1202     linesearch->rtol = rtol;
1203   }
1204 
1205   if (atol != (PetscReal)PETSC_DEFAULT) {
1206     PetscCheck(atol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Absolute tolerance %14.12e must be non-negative", (double)atol);
1207     linesearch->atol = atol;
1208   }
1209 
1210   if (ltol != (PetscReal)PETSC_DEFAULT) {
1211     PetscCheck(ltol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Lambda tolerance %14.12e must be non-negative", (double)ltol);
1212     linesearch->ltol = ltol;
1213   }
1214 
1215   if (max_it != PETSC_DEFAULT) {
1216     PetscCheck(max_it >= 0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of iterations %" PetscInt_FMT " must be non-negative", max_it);
1217     linesearch->max_its = max_it;
1218   }
1219   PetscFunctionReturn(PETSC_SUCCESS);
1220 }
1221 
1222 /*@
1223   SNESLineSearchGetDamping - Gets the line search damping parameter.
1224 
1225   Input Parameter:
1226 . linesearch - the line search context
1227 
1228   Output Parameter:
1229 . damping - The damping parameter
1230 
1231   Level: advanced
1232 
1233 .seealso: [](ch_snes), `SNES`, `SNESLineSearchGetStepTolerance()`, `SNESQN`
1234 @*/
1235 PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch linesearch, PetscReal *damping)
1236 {
1237   PetscFunctionBegin;
1238   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1239   PetscAssertPointer(damping, 2);
1240   *damping = linesearch->damping;
1241   PetscFunctionReturn(PETSC_SUCCESS);
1242 }
1243 
1244 /*@
1245   SNESLineSearchSetDamping - Sets the line search damping parameter.
1246 
1247   Input Parameters:
1248 + linesearch - the line search context
1249 - damping    - The damping parameter
1250 
1251   Options Database Key:
1252 . -snes_linesearch_damping <damping> - the damping value
1253 
1254   Level: intermediate
1255 
1256   Note:
1257   The `SNESLINESEARCHNONE` line search merely takes the update step scaled by the damping parameter.
1258   The use of the damping parameter in the `SNESLINESEARCHL2` and `SNESLINESEARCHCP` line searches is much more subtle;
1259   it is used as a starting point in calculating the secant step. However, the eventual
1260   step may be of greater length than the damping parameter.  In the `SNESLINESEARCHBT` line search it is
1261   used as the maximum possible step length, as the `SNESLINESEARCHBT` line search only backtracks.
1262 
1263 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetDamping()`
1264 @*/
1265 PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch linesearch, PetscReal damping)
1266 {
1267   PetscFunctionBegin;
1268   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1269   linesearch->damping = damping;
1270   PetscFunctionReturn(PETSC_SUCCESS);
1271 }
1272 
1273 /*@
1274   SNESLineSearchGetOrder - Gets the line search approximation order.
1275 
1276   Input Parameter:
1277 . linesearch - the line search context
1278 
1279   Output Parameter:
1280 . order - The order
1281 
1282   Level: intermediate
1283 
1284 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetOrder()`
1285 @*/
1286 PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch linesearch, PetscInt *order)
1287 {
1288   PetscFunctionBegin;
1289   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1290   PetscAssertPointer(order, 2);
1291   *order = linesearch->order;
1292   PetscFunctionReturn(PETSC_SUCCESS);
1293 }
1294 
1295 /*@
1296   SNESLineSearchSetOrder - Sets the maximum order of the polynomial fit used in the line search
1297 
1298   Input Parameters:
1299 + linesearch - the line search context
1300 - order      - The order
1301 
1302   Level: intermediate
1303 
1304   Values for `order`\:
1305 +  1 or `SNES_LINESEARCH_ORDER_LINEAR` - linear order
1306 .  2 or `SNES_LINESEARCH_ORDER_QUADRATIC` - quadratic order
1307 -  3 or `SNES_LINESEARCH_ORDER_CUBIC` - cubic order
1308 
1309   Options Database Key:
1310 . -snes_linesearch_order <order> - 1, 2, 3.  Most types only support certain orders (`SNESLINESEARCHBT` supports 2 or 3)
1311 
1312   Note:
1313   These orders are supported by `SNESLINESEARCHBT` and `SNESLINESEARCHCP`
1314 
1315 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetOrder()`, `SNESLineSearchSetDamping()`
1316 @*/
1317 PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch linesearch, PetscInt order)
1318 {
1319   PetscFunctionBegin;
1320   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1321   linesearch->order = order;
1322   PetscFunctionReturn(PETSC_SUCCESS);
1323 }
1324 
1325 /*@
1326   SNESLineSearchGetNorms - Gets the norms for the current solution `X`, the current update `Y`, and the current function value `F`.
1327 
1328   Not Collective
1329 
1330   Input Parameter:
1331 . linesearch - the line search context
1332 
1333   Output Parameters:
1334 + xnorm - The norm of the current solution
1335 . fnorm - The norm of the current function, this is the `norm(function(X))` where `X` is the current solution.
1336 - ynorm - The norm of the current update (after scaling by the linesearch computed lambda)
1337 
1338   Level: developer
1339 
1340   Notes:
1341   Some values may not be up-to-date at particular points in the code.
1342 
1343   This, in combination with `SNESLineSearchSetNorms()`, allow the line search and the `SNESSolve_XXX()` to share
1344   computed values.
1345 
1346 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetNorms()`, `SNESLineSearchGetVecs()`
1347 @*/
1348 PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal *xnorm, PetscReal *fnorm, PetscReal *ynorm)
1349 {
1350   PetscFunctionBegin;
1351   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1352   if (xnorm) *xnorm = linesearch->xnorm;
1353   if (fnorm) *fnorm = linesearch->fnorm;
1354   if (ynorm) *ynorm = linesearch->ynorm;
1355   PetscFunctionReturn(PETSC_SUCCESS);
1356 }
1357 
1358 /*@
1359   SNESLineSearchSetNorms - Sets the computed norms for the current solution `X`, the current update `Y`, and the current function value `F`.
1360 
1361   Collective
1362 
1363   Input Parameters:
1364 + linesearch - the line search context
1365 . xnorm      - The norm of the current solution
1366 . fnorm      - The norm of the current function, this is the `norm(function(X))` where `X` is the current solution
1367 - ynorm      - The norm of the current update (after scaling by the linesearch computed lambda)
1368 
1369   Level: developer
1370 
1371   Note:
1372   This is called by the line search routines to store the values they have just computed
1373 
1374 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()`
1375 @*/
1376 PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
1377 {
1378   PetscFunctionBegin;
1379   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1380   linesearch->xnorm = xnorm;
1381   linesearch->fnorm = fnorm;
1382   linesearch->ynorm = ynorm;
1383   PetscFunctionReturn(PETSC_SUCCESS);
1384 }
1385 
1386 /*@
1387   SNESLineSearchComputeNorms - Explicitly computes the norms of the current solution `X`, the current update `Y`, and the current function value `F`.
1388 
1389   Input Parameter:
1390 . linesearch - the line search context
1391 
1392   Options Database Key:
1393 . -snes_linesearch_norms - turn norm computation on or off
1394 
1395   Level: intermediate
1396 
1397   Developer Note:
1398   The options database key is misnamed. It should be -snes_linesearch_compute_norms
1399 
1400 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms`, `SNESLineSearchSetNorms()`, `SNESLineSearchSetComputeNorms()`
1401 @*/
1402 PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
1403 {
1404   SNES snes;
1405 
1406   PetscFunctionBegin;
1407   if (linesearch->norms) {
1408     if (linesearch->ops->vinorm) {
1409       PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
1410       PetscCall(VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
1411       PetscCall(VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm));
1412       PetscCall((*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm));
1413     } else {
1414       PetscCall(VecNormBegin(linesearch->vec_func, NORM_2, &linesearch->fnorm));
1415       PetscCall(VecNormBegin(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
1416       PetscCall(VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm));
1417       PetscCall(VecNormEnd(linesearch->vec_func, NORM_2, &linesearch->fnorm));
1418       PetscCall(VecNormEnd(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
1419       PetscCall(VecNormEnd(linesearch->vec_update, NORM_2, &linesearch->ynorm));
1420     }
1421   }
1422   PetscFunctionReturn(PETSC_SUCCESS);
1423 }
1424 
1425 /*@
1426   SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search.
1427 
1428   Input Parameters:
1429 + linesearch - the line search context
1430 - flg        - indicates whether or not to compute norms
1431 
1432   Options Database Key:
1433 . -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic (none) `SNESLINESEARCHBASIC` line search
1434 
1435   Level: intermediate
1436 
1437   Note:
1438   This is most relevant to the `SNESLINESEARCHBASIC` (or equivalently `SNESLINESEARCHNONE`) line search type since most line searches have a stopping criteria involving the norm.
1439 
1440   Developer Note:
1441   The options database key is misnamed. It should be -snes_linesearch_compute_norms
1442 
1443 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetNorms()`, `SNESLineSearchComputeNorms()`, `SNESLINESEARCHBASIC`
1444 @*/
1445 PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg)
1446 {
1447   PetscFunctionBegin;
1448   linesearch->norms = flg;
1449   PetscFunctionReturn(PETSC_SUCCESS);
1450 }
1451 
1452 /*@
1453   SNESLineSearchGetVecs - Gets the vectors from the `SNESLineSearch` context
1454 
1455   Not Collective but the vectors are parallel
1456 
1457   Input Parameter:
1458 . linesearch - the line search context
1459 
1460   Output Parameters:
1461 + X - Solution vector
1462 . F - Function vector
1463 . Y - Search direction vector
1464 . W - Solution work vector
1465 - G - Function work vector
1466 
1467   Level: advanced
1468 
1469   Notes:
1470   At the beginning of a line search application, `X` should contain a
1471   solution and the vector `F` the function computed at `X`.  At the end of the
1472   line search application, `X` should contain the new solution, and `F` the
1473   function evaluated at the new solution.
1474 
1475   These vectors are owned by the `SNESLineSearch` and should not be destroyed by the caller
1476 
1477 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()`
1478 @*/
1479 PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch, Vec *X, Vec *F, Vec *Y, Vec *W, Vec *G)
1480 {
1481   PetscFunctionBegin;
1482   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1483   if (X) {
1484     PetscAssertPointer(X, 2);
1485     *X = linesearch->vec_sol;
1486   }
1487   if (F) {
1488     PetscAssertPointer(F, 3);
1489     *F = linesearch->vec_func;
1490   }
1491   if (Y) {
1492     PetscAssertPointer(Y, 4);
1493     *Y = linesearch->vec_update;
1494   }
1495   if (W) {
1496     PetscAssertPointer(W, 5);
1497     *W = linesearch->vec_sol_new;
1498   }
1499   if (G) {
1500     PetscAssertPointer(G, 6);
1501     *G = linesearch->vec_func_new;
1502   }
1503   PetscFunctionReturn(PETSC_SUCCESS);
1504 }
1505 
1506 /*@
1507   SNESLineSearchSetVecs - Sets the vectors on the `SNESLineSearch` context
1508 
1509   Logically Collective
1510 
1511   Input Parameters:
1512 + linesearch - the line search context
1513 . X          - Solution vector
1514 . F          - Function vector
1515 . Y          - Search direction vector
1516 . W          - Solution work vector
1517 - G          - Function work vector
1518 
1519   Level: developer
1520 
1521 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetNorms()`, `SNESLineSearchGetVecs()`
1522 @*/
1523 PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch, Vec X, Vec F, Vec Y, Vec W, Vec G)
1524 {
1525   PetscFunctionBegin;
1526   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1527   if (X) {
1528     PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
1529     linesearch->vec_sol = X;
1530   }
1531   if (F) {
1532     PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
1533     linesearch->vec_func = F;
1534   }
1535   if (Y) {
1536     PetscValidHeaderSpecific(Y, VEC_CLASSID, 4);
1537     linesearch->vec_update = Y;
1538   }
1539   if (W) {
1540     PetscValidHeaderSpecific(W, VEC_CLASSID, 5);
1541     linesearch->vec_sol_new = W;
1542   }
1543   if (G) {
1544     PetscValidHeaderSpecific(G, VEC_CLASSID, 6);
1545     linesearch->vec_func_new = G;
1546   }
1547   PetscFunctionReturn(PETSC_SUCCESS);
1548 }
1549 
1550 /*@
1551   SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1552   `SNESLineSearch` options in the database.
1553 
1554   Logically Collective
1555 
1556   Input Parameters:
1557 + linesearch - the `SNESLineSearch` context
1558 - prefix     - the prefix to prepend to all option names
1559 
1560   Level: advanced
1561 
1562   Note:
1563   A hyphen (-) must NOT be given at the beginning of the prefix name.
1564   The first character of all runtime options is AUTOMATICALLY the hyphen.
1565 
1566 .seealso: [](ch_snes), `SNES`, `SNESLineSearch()`, `SNESLineSearchSetFromOptions()`, `SNESGetOptionsPrefix()`
1567 @*/
1568 PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch, const char prefix[])
1569 {
1570   PetscFunctionBegin;
1571   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1572   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)linesearch, prefix));
1573   PetscFunctionReturn(PETSC_SUCCESS);
1574 }
1575 
1576 /*@
1577   SNESLineSearchGetOptionsPrefix - Gets the prefix used for searching for all
1578   SNESLineSearch options in the database.
1579 
1580   Not Collective
1581 
1582   Input Parameter:
1583 . linesearch - the `SNESLineSearch` context
1584 
1585   Output Parameter:
1586 . prefix - pointer to the prefix string used
1587 
1588   Level: advanced
1589 
1590 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESAppendOptionsPrefix()`
1591 @*/
1592 PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch, const char *prefix[])
1593 {
1594   PetscFunctionBegin;
1595   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1596   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)linesearch, prefix));
1597   PetscFunctionReturn(PETSC_SUCCESS);
1598 }
1599 
1600 /*@C
1601   SNESLineSearchSetWorkVecs - Sets work vectors for the line search.
1602 
1603   Input Parameters:
1604 + linesearch - the `SNESLineSearch` context
1605 - nwork      - the number of work vectors
1606 
1607   Level: developer
1608 
1609   Developer Note:
1610   This is called from within the set up routines for each of the line search types `SNESLineSearchType`
1611 
1612 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESSetWorkVecs()`
1613 @*/
1614 PetscErrorCode SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork)
1615 {
1616   PetscFunctionBegin;
1617   PetscCheck(linesearch->vec_sol, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
1618   PetscCall(VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work));
1619   PetscFunctionReturn(PETSC_SUCCESS);
1620 }
1621 
1622 /*@
1623   SNESLineSearchGetReason - Gets the success/failure status of the last line search application
1624 
1625   Input Parameter:
1626 . linesearch - the line search context
1627 
1628   Output Parameter:
1629 . result - The success or failure status
1630 
1631   Level: developer
1632 
1633   Note:
1634   This is typically called after `SNESLineSearchApply()` in order to determine if the line search failed
1635   (and set into the `SNES` convergence accordingly).
1636 
1637 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetReason()`, `SNESLineSearchReason`
1638 @*/
1639 PetscErrorCode SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result)
1640 {
1641   PetscFunctionBegin;
1642   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1643   PetscAssertPointer(result, 2);
1644   *result = linesearch->result;
1645   PetscFunctionReturn(PETSC_SUCCESS);
1646 }
1647 
1648 /*@
1649   SNESLineSearchSetReason - Sets the success/failure status of the line search application
1650 
1651   Logically Collective; No Fortran Support
1652 
1653   Input Parameters:
1654 + linesearch - the line search context
1655 - result     - The success or failure status
1656 
1657   Level: developer
1658 
1659   Note:
1660   This is typically called in a `SNESLineSearchType` implementation of `SNESLineSearchApply()` or a `SNESLINESEARCHSHELL` implementation to set
1661   the success or failure of the line search method.
1662 
1663 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchReason`, `SNESLineSearchGetSResult()`
1664 @*/
1665 PetscErrorCode SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result)
1666 {
1667   PetscFunctionBegin;
1668   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1669   linesearch->result = result;
1670   PetscFunctionReturn(PETSC_SUCCESS);
1671 }
1672 
1673 // PetscClangLinter pragma disable: -fdoc-param-list-func-parameter-documentation
1674 /*@C
1675   SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.
1676 
1677   Logically Collective
1678 
1679   Input Parameters:
1680 + linesearch   - the linesearch object
1681 . projectfunc  - function for projecting the function to the bounds, see `SNESLineSearchVIProjectFn` for calling sequence
1682 . normfunc     - function for computing the norm of an active set, see `SNESLineSearchVINormFn` for calling sequence
1683 - dirderivfunc - function for computing the directional derivative of an active set, see `SNESLineSearchVIDirDerivFn` for calling sequence
1684 
1685   Level: advanced
1686 
1687   Notes:
1688   The VI solvers require projection of the solution to the feasible set.  `projectfunc` should implement this.
1689 
1690   The VI solvers require special evaluation of the function norm such that the norm is only calculated
1691   on the inactive set.  This should be implemented by `normfunc`.
1692 
1693   The VI solvers further require special evaluation of the directional derivative (when assuming that there exists some $G(x)$
1694   for which the `SNESFunctionFn` $F(x) = grad G(x)$) such that it is only calculated on the inactive set.
1695   This should be implemented by `dirderivfunc`.
1696 
1697 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetVIFunctions()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()`,
1698           `SNESLineSearchVIProjectFn`, `SNESLineSearchVINormFn`, `SNESLineSearchVIDirDerivFn`
1699 @*/
1700 PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFn *projectfunc, SNESLineSearchVINormFn *normfunc, SNESLineSearchVIDirDerivFn *dirderivfunc)
1701 {
1702   PetscFunctionBegin;
1703   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1704   if (projectfunc) linesearch->ops->viproject = projectfunc;
1705   if (normfunc) linesearch->ops->vinorm = normfunc;
1706   if (dirderivfunc) linesearch->ops->vidirderiv = dirderivfunc;
1707   PetscFunctionReturn(PETSC_SUCCESS);
1708 }
1709 
1710 /*@C
1711   SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.
1712 
1713   Not Collective
1714 
1715   Input Parameter:
1716 . linesearch - the line search context, obtain with `SNESGetLineSearch()`
1717 
1718   Output Parameters:
1719 + projectfunc  - function for projecting the function to the bounds, see `SNESLineSearchVIProjectFn` for calling sequence
1720 . normfunc     - function for computing the norm of an active set, see `SNESLineSearchVINormFn ` for calling sequence
1721 - dirderivfunc - function for computing the directional derivative of an active set, see `SNESLineSearchVIDirDerivFn` for calling sequence
1722 
1723   Level: advanced
1724 
1725 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetVIFunctions()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()`,
1726           `SNESLineSearchVIProjectFn`, `SNESLineSearchVINormFn`
1727 @*/
1728 PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFn **projectfunc, SNESLineSearchVINormFn **normfunc, SNESLineSearchVIDirDerivFn **dirderivfunc)
1729 {
1730   PetscFunctionBegin;
1731   if (projectfunc) *projectfunc = linesearch->ops->viproject;
1732   if (normfunc) *normfunc = linesearch->ops->vinorm;
1733   if (dirderivfunc) *dirderivfunc = linesearch->ops->vidirderiv;
1734   PetscFunctionReturn(PETSC_SUCCESS);
1735 }
1736 
1737 /*@C
1738   SNESLineSearchRegister - register a line search type `SNESLineSearchType`
1739 
1740   Logically Collective, No Fortran Support
1741 
1742   Input Parameters:
1743 + sname    - name of the `SNESLineSearchType()`
1744 - function - the creation function for that type
1745 
1746   Calling sequence of `function`:
1747 . ls - the line search context
1748 
1749   Level: advanced
1750 
1751 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchSetType()`
1752 @*/
1753 PetscErrorCode SNESLineSearchRegister(const char sname[], PetscErrorCode (*function)(SNESLineSearch ls))
1754 {
1755   PetscFunctionBegin;
1756   PetscCall(SNESInitializePackage());
1757   PetscCall(PetscFunctionListAdd(&SNESLineSearchList, sname, function));
1758   PetscFunctionReturn(PETSC_SUCCESS);
1759 }
1760