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