xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision b9c875b8ccb0466cfbdbcc08c79b487b3a7d395f)
1 #include <petsc/private/linesearchimpl.h> /*I "petscsnes.h" I*/
2 
3 PetscBool         SNESLineSearchRegisterAllCalled = PETSC_FALSE;
4 PetscFunctionList SNESLineSearchList              = NULL;
5 
6 PetscClassId  SNESLINESEARCH_CLASSID;
7 PetscLogEvent SNESLINESEARCH_Apply;
8 
9 /*@
10   SNESLineSearchMonitorCancel - Clears all the monitor functions for a `SNESLineSearch` object.
11 
12   Logically Collective
13 
14   Input Parameter:
15 . ls - the `SNESLineSearch` context
16 
17   Options Database Key:
18 . -snes_linesearch_monitor_cancel - cancels all monitors that have been hardwired
19     into a code by calls to `SNESLineSearchMonitorSet()`, but does not cancel those
20     set via the options database
21 
22   Level: advanced
23 
24   Notes:
25   There is no way to clear one specific monitor from a `SNESLineSearch` object.
26 
27   This does not clear the monitor set with `SNESLineSearchSetDefaultMonitor()` use `SNESLineSearchSetDefaultMonitor`(`ls`,`NULL`) to cancel it
28   that one.
29 
30 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorDefault()`, `SNESLineSearchMonitorSet()`
31 @*/
32 PetscErrorCode SNESLineSearchMonitorCancel(SNESLineSearch ls)
33 {
34   PetscInt i;
35 
36   PetscFunctionBegin;
37   PetscValidHeaderSpecific(ls, SNESLINESEARCH_CLASSID, 1);
38   for (i = 0; i < ls->numbermonitors; i++) {
39     if (ls->monitordestroy[i]) PetscCall((*ls->monitordestroy[i])(&ls->monitorcontext[i]));
40   }
41   ls->numbermonitors = 0;
42   PetscFunctionReturn(PETSC_SUCCESS);
43 }
44 
45 /*@
46   SNESLineSearchMonitor - runs the user provided monitor routines, if they exist
47 
48   Collective
49 
50   Input Parameter:
51 . ls - the linesearch object
52 
53   Level: developer
54 
55   Note:
56   This routine is called by the `SNESLineSearch` implementations.
57   It does not typically need to be called by the user.
58 
59 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorSet()`
60 @*/
61 PetscErrorCode SNESLineSearchMonitor(SNESLineSearch ls)
62 {
63   PetscInt i, n = ls->numbermonitors;
64 
65   PetscFunctionBegin;
66   for (i = 0; i < n; i++) PetscCall((*ls->monitorftns[i])(ls, ls->monitorcontext[i]));
67   PetscFunctionReturn(PETSC_SUCCESS);
68 }
69 
70 /*@C
71   SNESLineSearchMonitorSet - Sets an ADDITIONAL function that is to be used at every
72   iteration of the nonlinear solver to display the iteration's
73   progress.
74 
75   Logically Collective
76 
77   Input Parameters:
78 + ls             - the `SNESLineSearch` context
79 . f              - the monitor function
80 . mctx           - [optional] user-defined context for private data for the monitor routine (use `NULL` if no context is desired)
81 - monitordestroy - [optional] routine that frees monitor context (may be `NULL`)
82 
83   Calling sequence of `f`:
84 + ls   - the `SNESLineSearch` context
85 - mctx - [optional] user-defined context for private data for the monitor routine
86 
87   Calling sequence of `monitordestroy`:
88 . mctx - [optional] user-defined context for private data for the monitor routine
89 
90   Level: intermediate
91 
92   Note:
93   Several different monitoring routines may be set by calling
94   `SNESLineSearchMonitorSet()` multiple times; all will be called in the
95   order in which they were set.
96 
97   Fortran Note:
98   Only a single monitor function can be set for each `SNESLineSearch` object
99 
100 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorDefault()`, `SNESLineSearchMonitorCancel()`
101 @*/
102 PetscErrorCode SNESLineSearchMonitorSet(SNESLineSearch ls, PetscErrorCode (*f)(SNESLineSearch ls, void *mctx), void *mctx, PetscErrorCode (*monitordestroy)(void **mctx))
103 {
104   PetscInt  i;
105   PetscBool identical;
106 
107   PetscFunctionBegin;
108   PetscValidHeaderSpecific(ls, SNESLINESEARCH_CLASSID, 1);
109   for (i = 0; i < ls->numbermonitors; i++) {
110     PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))f, mctx, monitordestroy, (PetscErrorCode(*)(void))ls->monitorftns[i], ls->monitorcontext[i], ls->monitordestroy[i], &identical));
111     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
112   }
113   PetscCheck(ls->numbermonitors < MAXSNESLSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many monitors set");
114   ls->monitorftns[ls->numbermonitors]      = f;
115   ls->monitordestroy[ls->numbermonitors]   = monitordestroy;
116   ls->monitorcontext[ls->numbermonitors++] = (void *)mctx;
117   PetscFunctionReturn(PETSC_SUCCESS);
118 }
119 
120 /*@C
121   SNESLineSearchMonitorSolutionUpdate - Monitors each update of the function value the linesearch tries
122 
123   Collective
124 
125   Input Parameters:
126 + ls - the `SNESLineSearch` object
127 - vf - the context for the monitor, in this case it is an `PetscViewerAndFormat`
128 
129   Options Database Key:
130 . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine
131 
132   Level: developer
133 
134   This is not normally called directly but is passed to `SNESLineSearchMonitorSet()`
135 
136 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorSet()`, `SNESMonitorSolution()`
137 @*/
138 PetscErrorCode SNESLineSearchMonitorSolutionUpdate(SNESLineSearch ls, PetscViewerAndFormat *vf)
139 {
140   PetscViewer viewer = vf->viewer;
141   Vec         Y, W, G;
142 
143   PetscFunctionBegin;
144   PetscCall(SNESLineSearchGetVecs(ls, NULL, NULL, &Y, &W, &G));
145   PetscCall(PetscViewerPushFormat(viewer, vf->format));
146   PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted update to solution \n"));
147   PetscCall(VecView(Y, viewer));
148   PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted new solution \n"));
149   PetscCall(VecView(W, viewer));
150   PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted updated function value\n"));
151   PetscCall(VecView(G, viewer));
152   PetscCall(PetscViewerPopFormat(viewer));
153   PetscFunctionReturn(PETSC_SUCCESS);
154 }
155 
156 /*@
157   SNESLineSearchCreate - Creates a `SNESLineSearch` context.
158 
159   Logically Collective
160 
161   Input Parameter:
162 . comm - MPI communicator for the line search (typically from the associated `SNES` context).
163 
164   Output Parameter:
165 . outlinesearch - the new line search context
166 
167   Level: developer
168 
169   Note:
170   The preferred calling sequence is to use `SNESGetLineSearch()` to acquire the `SNESLineSearch` instance
171   already associated with the `SNES`.
172 
173 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `LineSearchDestroy()`, `SNESGetLineSearch()`
174 @*/
175 PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch)
176 {
177   SNESLineSearch linesearch;
178 
179   PetscFunctionBegin;
180   PetscAssertPointer(outlinesearch, 2);
181   PetscCall(SNESInitializePackage());
182   *outlinesearch = NULL;
183 
184   PetscCall(PetscHeaderCreate(linesearch, SNESLINESEARCH_CLASSID, "SNESLineSearch", "Linesearch", "SNESLineSearch", comm, SNESLineSearchDestroy, SNESLineSearchView));
185 
186   linesearch->vec_sol_new  = NULL;
187   linesearch->vec_func_new = NULL;
188   linesearch->vec_sol      = NULL;
189   linesearch->vec_func     = NULL;
190   linesearch->vec_update   = NULL;
191 
192   linesearch->lambda       = 1.0;
193   linesearch->fnorm        = 1.0;
194   linesearch->ynorm        = 1.0;
195   linesearch->xnorm        = 1.0;
196   linesearch->result       = SNES_LINESEARCH_SUCCEEDED;
197   linesearch->norms        = PETSC_TRUE;
198   linesearch->keeplambda   = PETSC_FALSE;
199   linesearch->damping      = 1.0;
200   linesearch->maxstep      = 1e8;
201   linesearch->steptol      = 1e-12;
202   linesearch->rtol         = 1e-8;
203   linesearch->atol         = 1e-15;
204   linesearch->ltol         = 1e-8;
205   linesearch->precheckctx  = NULL;
206   linesearch->postcheckctx = NULL;
207   linesearch->max_its      = 1;
208   linesearch->setupcalled  = PETSC_FALSE;
209   linesearch->monitor      = NULL;
210   *outlinesearch           = linesearch;
211   PetscFunctionReturn(PETSC_SUCCESS);
212 }
213 
214 /*@
215   SNESLineSearchSetUp - Prepares the line search for being applied by allocating
216   any required vectors.
217 
218   Collective
219 
220   Input Parameter:
221 . linesearch - The `SNESLineSearch` instance.
222 
223   Level: advanced
224 
225   Note:
226   For most cases, this needn't be called by users or outside of `SNESLineSearchApply()`.
227   The only current case where this is called outside of this is for the VI
228   solvers, which modify the solution and work vectors before the first call
229   of `SNESLineSearchApply()`, requiring the `SNESLineSearch` work vectors to be
230   allocated upfront.
231 
232 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchReset()`
233 @*/
234 PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch)
235 {
236   PetscFunctionBegin;
237   if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC));
238   if (!linesearch->setupcalled) {
239     if (!linesearch->vec_sol_new) PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new));
240     if (!linesearch->vec_func_new) PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new));
241     PetscTryTypeMethod(linesearch, setup);
242     if (!linesearch->ops->snesfunc) PetscCall(SNESLineSearchSetFunction(linesearch, SNESComputeFunction));
243     linesearch->lambda      = linesearch->damping;
244     linesearch->setupcalled = PETSC_TRUE;
245   }
246   PetscFunctionReturn(PETSC_SUCCESS);
247 }
248 
249 /*@
250   SNESLineSearchReset - Undoes the `SNESLineSearchSetUp()` and deletes any `Vec`s or `Mat`s allocated by the line search.
251 
252   Collective
253 
254   Input Parameter:
255 . linesearch - The `SNESLineSearch` instance.
256 
257   Level: developer
258 
259   Note:
260   Usually only called by `SNESReset()`
261 
262 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchSetUp()`
263 @*/
264 PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch)
265 {
266   PetscFunctionBegin;
267   PetscTryTypeMethod(linesearch, reset);
268 
269   PetscCall(VecDestroy(&linesearch->vec_sol_new));
270   PetscCall(VecDestroy(&linesearch->vec_func_new));
271 
272   PetscCall(VecDestroyVecs(linesearch->nwork, &linesearch->work));
273 
274   linesearch->nwork       = 0;
275   linesearch->setupcalled = PETSC_FALSE;
276   PetscFunctionReturn(PETSC_SUCCESS);
277 }
278 
279 /*@C
280   SNESLineSearchSetFunction - Sets the function evaluation used by the `SNES` line search
281   `
282 
283   Input Parameters:
284 + linesearch - the `SNESLineSearch` context
285 - func       - function evaluation routine, this is usually the function provided with `SNESSetFunction()`
286 
287   Calling sequence of `func`:
288 + snes - the `SNES` with which the `SNESLineSearch` context is associated with
289 . x    - the input vector
290 - f    - the computed value of the function
291 
292   Level: developer
293 
294   Note:
295   By default the `SNESLineSearch` uses the function provided by `SNESSetFunction()` so this is rarely needed
296 
297 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESSetFunction()`
298 @*/
299 PetscErrorCode SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES snes, Vec x, Vec f))
300 {
301   PetscFunctionBegin;
302   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
303   linesearch->ops->snesfunc = func;
304   PetscFunctionReturn(PETSC_SUCCESS);
305 }
306 
307 /*@C
308   SNESLineSearchSetPreCheck - Sets a function that is called after the initial search direction has been computed but
309   before the line search routine has been applied. Allows adjusting the result of (usually a linear solve) that
310   determined the search direction.
311 
312   Logically Collective
313 
314   Input Parameters:
315 + linesearch - the `SNESLineSearch` context
316 . func       - [optional] function evaluation routine
317 - ctx        - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
318 
319   Calling sequence of `func`:
320 + ls        - the `SNESLineSearch` context
321 . x         - the current solution
322 . d         - the current search direction
323 . changed_d - indicates if the search direction has been changed
324 - ctx       - the context passed to `SNESLineSearchSetPreCheck()`
325 
326   Level: intermediate
327 
328   Note:
329   Use `SNESLineSearchSetPostCheck()` to change the step after the line search is complete.
330 
331   Use `SNESVISetVariableBounds()` and `SNESVISetComputeVariableBounds()` to cause `SNES` to automatically control the ranges of variables allowed.
332 
333 .seealso: [](ch_snes), `SNES`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()`,
334           `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`
335 
336 @*/
337 PetscErrorCode SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch ls, Vec x, Vec d, PetscBool *changed_d, void *ctx), void *ctx)
338 {
339   PetscFunctionBegin;
340   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
341   if (func) linesearch->ops->precheck = func;
342   if (ctx) linesearch->precheckctx = ctx;
343   PetscFunctionReturn(PETSC_SUCCESS);
344 }
345 
346 /*@C
347   SNESLineSearchGetPreCheck - Gets the pre-check function for the line search routine.
348 
349   Input Parameter:
350 . linesearch - the `SNESLineSearch` context
351 
352   Output Parameters:
353 + func - [optional] function evaluation routine,  for calling sequence see `SNESLineSearchSetPreCheck()`
354 - ctx  - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
355 
356   Level: intermediate
357 
358 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchSetPostCheck()`
359 @*/
360 PetscErrorCode SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch, Vec, Vec, PetscBool *, void *), void **ctx)
361 {
362   PetscFunctionBegin;
363   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
364   if (func) *func = linesearch->ops->precheck;
365   if (ctx) *ctx = linesearch->precheckctx;
366   PetscFunctionReturn(PETSC_SUCCESS);
367 }
368 
369 /*@C
370   SNESLineSearchSetPostCheck - Sets a user function that is called after the line search has been applied to determine the step
371   direction and length. Allows the user a chance to change or override the decision of the line search routine
372 
373   Logically Collective
374 
375   Input Parameters:
376 + linesearch - the `SNESLineSearch` context
377 . func       - [optional] function evaluation routine
378 - ctx        - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
379 
380   Calling sequence of `func`:
381 + ls        - the `SNESLineSearch` context
382 . x         - the current solution
383 . d         - the current search direction
384 . w         - $ w = x + lambda*d $ for some lambda
385 . changed_d - indicates if the search direction `d` has been changed
386 . changed_w - indicates `w` has been changed
387 - ctx       - the context passed to `SNESLineSearchSetPreCheck()`
388 
389   Level: intermediate
390 
391   Notes:
392   Use `SNESLineSearchSetPreCheck()` to change the step before the line search is completed.
393   The calling sequence of the callback does not contain the current scaling factor. To access the value, use `SNESLineSearchGetLambda()`.
394 
395   Use `SNESVISetVariableBounds()` and `SNESVISetComputeVariableBounds()` to cause `SNES` to automatically control the ranges of variables allowed.
396 
397 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchGetPostCheck()`,
398           `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`
399 @*/
400 PetscErrorCode SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch ls, Vec x, Vec d, Vec w, PetscBool *changed_d, PetscBool *changed_w, void *ctx), void *ctx)
401 {
402   PetscFunctionBegin;
403   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
404   if (func) linesearch->ops->postcheck = func;
405   if (ctx) linesearch->postcheckctx = ctx;
406   PetscFunctionReturn(PETSC_SUCCESS);
407 }
408 
409 /*@C
410   SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine.
411 
412   Input Parameter:
413 . linesearch - the `SNESLineSearch` context
414 
415   Output Parameters:
416 + func - [optional] function evaluation routine, see for the calling sequence `SNESLineSearchSetPostCheck()`
417 - ctx  - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
418 
419   Level: intermediate
420 
421 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`
422 @*/
423 PetscErrorCode SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx)
424 {
425   PetscFunctionBegin;
426   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
427   if (func) *func = linesearch->ops->postcheck;
428   if (ctx) *ctx = linesearch->postcheckctx;
429   PetscFunctionReturn(PETSC_SUCCESS);
430 }
431 
432 /*@
433   SNESLineSearchPreCheck - Prepares the line search for being applied.
434 
435   Logically Collective
436 
437   Input Parameters:
438 + linesearch - The linesearch instance.
439 . X          - The current solution
440 - Y          - The step direction
441 
442   Output Parameter:
443 . changed - Indicator that the precheck routine has changed `Y`
444 
445   Level: advanced
446 
447   Note:
448   This calls any function provided with `SNESLineSearchSetPreCheck()` and is called automatically inside the line search routines
449 
450   Developer Note:
451   The use of `PetscObjectGetState()` would eliminate the need for the `changed` argument to be provided
452 
453 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`,
454           `SNESLineSearchGetPostCheck()`
455 @*/
456 PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed)
457 {
458   PetscFunctionBegin;
459   *changed = PETSC_FALSE;
460   if (linesearch->ops->precheck) {
461     PetscUseTypeMethod(linesearch, precheck, X, Y, changed, linesearch->precheckctx);
462     PetscValidLogicalCollectiveBool(linesearch, *changed, 4);
463   }
464   PetscFunctionReturn(PETSC_SUCCESS);
465 }
466 
467 /*@
468   SNESLineSearchPostCheck - Hook to modify step direction or updated solution after a successful linesearch
469 
470   Logically Collective
471 
472   Input Parameters:
473 + linesearch - The line search context
474 . X          - The last solution
475 . Y          - The step direction
476 - W          - The updated solution, `W = X - lambda * Y` for some lambda
477 
478   Output Parameters:
479 + changed_Y - Indicator if the direction `Y` has been changed.
480 - changed_W - Indicator if the new candidate solution `W` has been changed.
481 
482   Level: developer
483 
484   Note:
485   This calls any function provided with `SNESLineSearchSetPostCheck()` and is called automatically inside the line search routines
486 
487   Developer Note:
488   The use of `PetscObjectGetState()` would eliminate the need for the `changed_Y` and `changed_W` arguments to be provided
489 
490 .seealso: [](ch_snes), `SNES`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPrecheck()`, `SNESLineSearchGetPrecheck()`
491 @*/
492 PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W)
493 {
494   PetscFunctionBegin;
495   *changed_Y = PETSC_FALSE;
496   *changed_W = PETSC_FALSE;
497   if (linesearch->ops->postcheck) {
498     PetscUseTypeMethod(linesearch, postcheck, X, Y, W, changed_Y, changed_W, linesearch->postcheckctx);
499     PetscValidLogicalCollectiveBool(linesearch, *changed_Y, 5);
500     PetscValidLogicalCollectiveBool(linesearch, *changed_W, 6);
501   }
502   PetscFunctionReturn(PETSC_SUCCESS);
503 }
504 
505 /*@C
506   SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration {cite}`hindmarsh1996time`
507 
508   Logically Collective
509 
510   Input Parameters:
511 + linesearch - the line search context
512 . X          - base state for this step
513 - ctx        - context for this function
514 
515   Input/Output Parameter:
516 . Y - correction, possibly modified
517 
518   Output Parameter:
519 . changed - flag indicating that `Y` was modified
520 
521   Options Database Keys:
522 + -snes_linesearch_precheck_picard       - activate this routine
523 - -snes_linesearch_precheck_picard_angle - angle
524 
525   Level: advanced
526 
527   Notes:
528   This function should be passed to `SNESLineSearchSetPreCheck()`
529 
530   The justification for this method involves the linear convergence of a Picard iteration
531   so the Picard linearization should be provided in place of the "Jacobian"  {cite}`hindmarsh1996time`. This correction
532   is generally not useful when using a Newton linearization.
533 
534   Developer Note:
535   The use of `PetscObjectGetState()` would eliminate the need for the `changed` argument to be provided
536 
537 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESSetPicard()`, `SNESGetLineSearch()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchSetPostCheck()`
538 @*/
539 PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed, void *ctx)
540 {
541   PetscReal   angle = *(PetscReal *)linesearch->precheckctx;
542   Vec         Ylast;
543   PetscScalar dot;
544   PetscInt    iter;
545   PetscReal   ynorm, ylastnorm, theta, angle_radians;
546   SNES        snes;
547 
548   PetscFunctionBegin;
549   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
550   PetscCall(PetscObjectQuery((PetscObject)snes, "SNESLineSearchPreCheckPicard_Ylast", (PetscObject *)&Ylast));
551   if (!Ylast) {
552     PetscCall(VecDuplicate(Y, &Ylast));
553     PetscCall(PetscObjectCompose((PetscObject)snes, "SNESLineSearchPreCheckPicard_Ylast", (PetscObject)Ylast));
554     PetscCall(PetscObjectDereference((PetscObject)Ylast));
555   }
556   PetscCall(SNESGetIterationNumber(snes, &iter));
557   if (iter < 2) {
558     PetscCall(VecCopy(Y, Ylast));
559     *changed = PETSC_FALSE;
560     PetscFunctionReturn(PETSC_SUCCESS);
561   }
562 
563   PetscCall(VecDot(Y, Ylast, &dot));
564   PetscCall(VecNorm(Y, NORM_2, &ynorm));
565   PetscCall(VecNorm(Ylast, NORM_2, &ylastnorm));
566   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
567   theta         = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm), -1.0, 1.0));
568   angle_radians = angle * PETSC_PI / 180.;
569   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
570     /* Modify the step Y */
571     PetscReal alpha, ydiffnorm;
572     PetscCall(VecAXPY(Ylast, -1.0, Y));
573     PetscCall(VecNorm(Ylast, NORM_2, &ydiffnorm));
574     alpha = (ydiffnorm > .001 * ylastnorm) ? ylastnorm / ydiffnorm : 1000.0;
575     PetscCall(VecCopy(Y, Ylast));
576     PetscCall(VecScale(Y, alpha));
577     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));
578     *changed = PETSC_TRUE;
579   } else {
580     PetscCall(PetscInfo(snes, "Angle %14.12e degrees exceeds threshold %14.12e, no correction applied\n", (double)(theta * 180. / PETSC_PI), (double)angle));
581     PetscCall(VecCopy(Y, Ylast));
582     *changed = PETSC_FALSE;
583   }
584   PetscFunctionReturn(PETSC_SUCCESS);
585 }
586 
587 /*@
588   SNESLineSearchApply - Computes the line-search update.
589 
590   Collective
591 
592   Input Parameter:
593 . linesearch - The line search context
594 
595   Input/Output Parameters:
596 + X     - The current solution, on output the new solution
597 . F     - The current function value, on output the new function value at the solution value `X`
598 . fnorm - The current norm of `F`, on output the new norm of `F`
599 - Y     - The current search direction, on output the direction determined by the linesearch, i.e. Xnew = Xold - lambda*Y
600 
601   Options Database Keys:
602 + -snes_linesearch_type                - basic (or equivalently none), bt, l2, cp, nleqerr, shell
603 . -snes_linesearch_monitor [:filename] - Print progress of line searches
604 . -snes_linesearch_damping             - The linesearch damping parameter, default is 1.0 (no damping)
605 . -snes_linesearch_norms               - Turn on/off the linesearch norms computation (SNESLineSearchSetComputeNorms())
606 . -snes_linesearch_keeplambda          - Keep the previous search length as the initial guess
607 - -snes_linesearch_max_it              - The number of iterations for iterative line searches
608 
609   Level: intermediate
610 
611   Notes:
612   This is typically called from within a `SNESSolve()` implementation in order to
613   help with convergence of the nonlinear method.  Various `SNES` types use line searches
614   in different ways, but the overarching theme is that a line search is used to determine
615   an optimal damping parameter of a step at each iteration of the method.  Each
616   application of the line search may invoke `SNESComputeFunction()` several times, and
617   therefore may be fairly expensive.
618 
619 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchGetLambda()`, `SNESLineSearchPreCheck()`, `SNESLineSearchPostCheck()`, `SNESSolve()`, `SNESComputeFunction()`, `SNESLineSearchSetComputeNorms()`,
620           `SNESLineSearchType`, `SNESLineSearchSetType()`
621 @*/
622 PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal *fnorm, Vec Y)
623 {
624   PetscFunctionBegin;
625   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
626   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
627   PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
628   PetscValidHeaderSpecific(Y, VEC_CLASSID, 5);
629 
630   linesearch->result = SNES_LINESEARCH_SUCCEEDED;
631 
632   linesearch->vec_sol    = X;
633   linesearch->vec_update = Y;
634   linesearch->vec_func   = F;
635 
636   PetscCall(SNESLineSearchSetUp(linesearch));
637 
638   if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */
639 
640   if (fnorm) linesearch->fnorm = *fnorm;
641   else PetscCall(VecNorm(F, NORM_2, &linesearch->fnorm));
642 
643   PetscCall(PetscLogEventBegin(SNESLINESEARCH_Apply, linesearch, X, F, Y));
644 
645   PetscUseTypeMethod(linesearch, apply);
646 
647   PetscCall(PetscLogEventEnd(SNESLINESEARCH_Apply, linesearch, X, F, Y));
648 
649   if (fnorm) *fnorm = linesearch->fnorm;
650   PetscFunctionReturn(PETSC_SUCCESS);
651 }
652 
653 /*@
654   SNESLineSearchDestroy - Destroys the line search instance.
655 
656   Collective
657 
658   Input Parameter:
659 . linesearch - The line search context
660 
661   Level: developer
662 
663   Note:
664   The line search in `SNES` is automatically called on `SNESDestroy()` so this call is rarely needed
665 
666 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchReset()`, `SNESDestroy()`
667 @*/
668 PetscErrorCode SNESLineSearchDestroy(SNESLineSearch *linesearch)
669 {
670   PetscFunctionBegin;
671   if (!*linesearch) PetscFunctionReturn(PETSC_SUCCESS);
672   PetscValidHeaderSpecific(*linesearch, SNESLINESEARCH_CLASSID, 1);
673   if (--((PetscObject)*linesearch)->refct > 0) {
674     *linesearch = NULL;
675     PetscFunctionReturn(PETSC_SUCCESS);
676   }
677   PetscCall(PetscObjectSAWsViewOff((PetscObject)*linesearch));
678   PetscCall(SNESLineSearchReset(*linesearch));
679   PetscTryTypeMethod(*linesearch, destroy);
680   PetscCall(PetscOptionsRestoreViewer(&(*linesearch)->monitor));
681   PetscCall(SNESLineSearchMonitorCancel(*linesearch));
682   PetscCall(PetscHeaderDestroy(linesearch));
683   PetscFunctionReturn(PETSC_SUCCESS);
684 }
685 
686 /*@
687   SNESLineSearchSetDefaultMonitor - Turns on/off printing useful information and debugging output about the line search.
688 
689   Logically Collective
690 
691   Input Parameters:
692 + linesearch - the linesearch object
693 - viewer     - an `PETSCVIEWERASCII` `PetscViewer` or `NULL` to turn off monitor
694 
695   Options Database Key:
696 . -snes_linesearch_monitor [:filename] - enables the monitor
697 
698   Level: intermediate
699 
700   Developer Notes:
701   This monitor is implemented differently than the other line search monitors that are set with
702   `SNESLineSearchMonitorSet()` since it is called in many locations of the line search routines to display aspects of the
703   line search that are not visible to the other monitors.
704 
705 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `PETSCVIEWERASCII`, `SNESGetLineSearch()`, `SNESLineSearchGetDefaultMonitor()`, `PetscViewer`, `SNESLineSearchSetMonitor()`,
706           `SNESLineSearchMonitorSetFromOptions()`
707 @*/
708 PetscErrorCode SNESLineSearchSetDefaultMonitor(SNESLineSearch linesearch, PetscViewer viewer)
709 {
710   PetscFunctionBegin;
711   PetscCall(PetscOptionsRestoreViewer(&linesearch->monitor));
712   linesearch->monitor = viewer;
713   PetscFunctionReturn(PETSC_SUCCESS);
714 }
715 
716 /*@
717   SNESLineSearchGetDefaultMonitor - Gets the `PetscViewer` instance for the default line search monitor that is turned on with `SNESLineSearchSetDefaultMonitor()`
718 
719   Logically Collective
720 
721   Input Parameter:
722 . linesearch - the line search context
723 
724   Output Parameter:
725 . monitor - monitor context
726 
727   Level: intermediate
728 
729 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchSetDefaultMonitor()`, `PetscViewer`
730 @*/
731 PetscErrorCode SNESLineSearchGetDefaultMonitor(SNESLineSearch linesearch, PetscViewer *monitor)
732 {
733   PetscFunctionBegin;
734   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
735   *monitor = linesearch->monitor;
736   PetscFunctionReturn(PETSC_SUCCESS);
737 }
738 
739 /*@C
740   SNESLineSearchMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated in the options database
741 
742   Collective
743 
744   Input Parameters:
745 + ls           - `SNESLineSearch` object to monitor
746 . name         - the monitor type
747 . help         - message indicating what monitoring is done
748 . manual       - manual page for the monitor
749 . monitor      - the monitor function
750 - monitorsetup - a function that is called once ONLY if the user selected this monitor that may set additional features of the `SNESLineSearch` or `PetscViewer`
751 
752   Calling sequence of `monitor`:
753 + ls - `SNESLineSearch` object being monitored
754 - vf - a `PetscViewerAndFormat` struct that provides the `PetscViewer` and `PetscViewerFormat` being used
755 
756   Calling sequence of `monitorsetup`:
757 + ls - `SNESLineSearch` object being monitored
758 - vf - a `PetscViewerAndFormat` struct that provides the `PetscViewer` and `PetscViewerFormat` being used
759 
760   Level: advanced
761 
762 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetMonitor()`, `PetscOptionsGetViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
763           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
764           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
765           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
766           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
767           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
768           `PetscOptionsFList()`, `PetscOptionsEList()`
769 @*/
770 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))
771 {
772   PetscViewer       viewer;
773   PetscViewerFormat format;
774   PetscBool         flg;
775 
776   PetscFunctionBegin;
777   PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)ls), ((PetscObject)ls)->options, ((PetscObject)ls)->prefix, name, &viewer, &format, &flg));
778   if (flg) {
779     PetscViewerAndFormat *vf;
780     PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf));
781     PetscCall(PetscOptionsRestoreViewer(&viewer));
782     if (monitorsetup) PetscCall((*monitorsetup)(ls, vf));
783     PetscCall(SNESLineSearchMonitorSet(ls, (PetscErrorCode(*)(SNESLineSearch, void *))monitor, vf, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy));
784   }
785   PetscFunctionReturn(PETSC_SUCCESS);
786 }
787 
788 /*@
789   SNESLineSearchSetFromOptions - Sets options for the line search
790 
791   Logically Collective
792 
793   Input Parameter:
794 . linesearch - a `SNESLineSearch` line search context
795 
796   Options Database Keys:
797 + -snes_linesearch_type <type>                                      - basic (or equivalently none), bt, l2, cp, nleqerr, shell
798 . -snes_linesearch_order <order>                                    - 1, 2, 3.  Most types only support certain orders (bt supports 2 or 3)
799 . -snes_linesearch_norms                                            - Turn on/off the linesearch norms for the basic linesearch typem (`SNESLineSearchSetComputeNorms()`)
800 . -snes_linesearch_minlambda                                        - The minimum step length
801 . -snes_linesearch_maxstep                                          - The maximum step size
802 . -snes_linesearch_rtol                                             - Relative tolerance for iterative line searches
803 . -snes_linesearch_atol                                             - Absolute tolerance for iterative line searches
804 . -snes_linesearch_ltol                                             - Change in lambda tolerance for iterative line searches
805 . -snes_linesearch_max_it                                           - The number of iterations for iterative line searches
806 . -snes_linesearch_monitor [:filename]                              - Print progress of line searches
807 . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine
808 . -snes_linesearch_damping                                          - The linesearch damping parameter
809 . -snes_linesearch_keeplambda                                       - Keep the previous search length as the initial guess.
810 . -snes_linesearch_precheck_picard                                  - Use precheck that speeds up convergence of picard method
811 - -snes_linesearch_precheck_picard_angle                            - Angle used in Picard precheck method
812 
813   Level: intermediate
814 
815 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchCreate()`, `SNESLineSearchSetOrder()`, `SNESLineSearchSetType()`, `SNESLineSearchSetTolerances()`, `SNESLineSearchSetDamping()`, `SNESLineSearchPreCheckPicard()`,
816           `SNESLineSearchType`, `SNESLineSearchSetComputeNorms()`
817 @*/
818 PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch)
819 {
820   const char *deft = SNESLINESEARCHBASIC;
821   char        type[256];
822   PetscBool   flg, set;
823   PetscViewer viewer;
824 
825   PetscFunctionBegin;
826   PetscCall(SNESLineSearchRegisterAll());
827 
828   PetscObjectOptionsBegin((PetscObject)linesearch);
829   if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name;
830   PetscCall(PetscOptionsFList("-snes_linesearch_type", "Linesearch type", "SNESLineSearchSetType", SNESLineSearchList, deft, type, 256, &flg));
831   if (flg) {
832     PetscCall(SNESLineSearchSetType(linesearch, type));
833   } else if (!((PetscObject)linesearch)->type_name) {
834     PetscCall(SNESLineSearchSetType(linesearch, deft));
835   }
836 
837   PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)linesearch), ((PetscObject)linesearch)->options, ((PetscObject)linesearch)->prefix, "-snes_linesearch_monitor", &viewer, NULL, &set));
838   if (set) {
839     PetscCall(SNESLineSearchSetDefaultMonitor(linesearch, viewer));
840     PetscCall(PetscOptionsRestoreViewer(&viewer));
841   }
842   PetscCall(SNESLineSearchMonitorSetFromOptions(linesearch, "-snes_linesearch_monitor_solution_update", "View correction at each iteration", "SNESLineSearchMonitorSolutionUpdate", SNESLineSearchMonitorSolutionUpdate, NULL));
843 
844   /* tolerances */
845   PetscCall(PetscOptionsReal("-snes_linesearch_minlambda", "Minimum step length", "SNESLineSearchSetTolerances", linesearch->steptol, &linesearch->steptol, NULL));
846   PetscCall(PetscOptionsReal("-snes_linesearch_maxstep", "Maximum step size", "SNESLineSearchSetTolerances", linesearch->maxstep, &linesearch->maxstep, NULL));
847   PetscCall(PetscOptionsReal("-snes_linesearch_rtol", "Relative tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->rtol, &linesearch->rtol, NULL));
848   PetscCall(PetscOptionsReal("-snes_linesearch_atol", "Absolute tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->atol, &linesearch->atol, NULL));
849   PetscCall(PetscOptionsReal("-snes_linesearch_ltol", "Change in lambda tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->ltol, &linesearch->ltol, NULL));
850   PetscCall(PetscOptionsInt("-snes_linesearch_max_it", "Maximum iterations for iterative line searches", "SNESLineSearchSetTolerances", linesearch->max_its, &linesearch->max_its, NULL));
851 
852   /* damping parameters */
853   PetscCall(PetscOptionsReal("-snes_linesearch_damping", "Line search damping and initial step guess", "SNESLineSearchSetDamping", linesearch->damping, &linesearch->damping, NULL));
854 
855   PetscCall(PetscOptionsBool("-snes_linesearch_keeplambda", "Use previous lambda as damping", "SNESLineSearchSetKeepLambda", linesearch->keeplambda, &linesearch->keeplambda, NULL));
856 
857   /* precheck */
858   PetscCall(PetscOptionsBool("-snes_linesearch_precheck_picard", "Use a correction that sometimes improves convergence of Picard iteration", "SNESLineSearchPreCheckPicard", flg, &flg, &set));
859   if (set) {
860     if (flg) {
861       linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */
862 
863       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));
864       PetscCall(SNESLineSearchSetPreCheck(linesearch, SNESLineSearchPreCheckPicard, &linesearch->precheck_picard_angle));
865     } else {
866       PetscCall(SNESLineSearchSetPreCheck(linesearch, NULL, NULL));
867     }
868   }
869   PetscCall(PetscOptionsInt("-snes_linesearch_order", "Order of approximation used in the line search", "SNESLineSearchSetOrder", linesearch->order, &linesearch->order, NULL));
870   PetscCall(PetscOptionsBool("-snes_linesearch_norms", "Compute final norms in line search", "SNESLineSearchSetComputeNorms", linesearch->norms, &linesearch->norms, NULL));
871 
872   PetscTryTypeMethod(linesearch, setfromoptions, PetscOptionsObject);
873 
874   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)linesearch, PetscOptionsObject));
875   PetscOptionsEnd();
876   PetscFunctionReturn(PETSC_SUCCESS);
877 }
878 
879 /*@
880   SNESLineSearchView - Prints useful information about the line search
881 
882   Logically Collective
883 
884   Input Parameters:
885 + linesearch - line search context
886 - viewer     - the `PetscViewer` to display the line search information to
887 
888   Level: intermediate
889 
890 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `PetscViewer`, `SNESLineSearchCreate()`
891 @*/
892 PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer)
893 {
894   PetscBool iascii;
895 
896   PetscFunctionBegin;
897   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
898   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch), &viewer));
899   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
900   PetscCheckSameComm(linesearch, 1, viewer, 2);
901 
902   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
903   if (iascii) {
904     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)linesearch, viewer));
905     PetscCall(PetscViewerASCIIPushTab(viewer));
906     PetscTryTypeMethod(linesearch, view, viewer);
907     PetscCall(PetscViewerASCIIPopTab(viewer));
908     PetscCall(PetscViewerASCIIPrintf(viewer, "  maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep, (double)linesearch->steptol));
909     PetscCall(PetscViewerASCIIPrintf(viewer, "  tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol, (double)linesearch->atol, (double)linesearch->ltol));
910     PetscCall(PetscViewerASCIIPrintf(viewer, "  maximum iterations=%" PetscInt_FMT "\n", linesearch->max_its));
911     if (linesearch->ops->precheck) {
912       if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) {
913         PetscCall(PetscViewerASCIIPrintf(viewer, "  using precheck step to speed up Picard convergence\n"));
914       } else {
915         PetscCall(PetscViewerASCIIPrintf(viewer, "  using user-defined precheck step\n"));
916       }
917     }
918     if (linesearch->ops->postcheck) PetscCall(PetscViewerASCIIPrintf(viewer, "  using user-defined postcheck step\n"));
919   }
920   PetscFunctionReturn(PETSC_SUCCESS);
921 }
922 
923 /*@C
924   SNESLineSearchGetType - Gets the `SNESLinesearchType` of a `SNESLineSearch`
925 
926   Logically Collective
927 
928   Input Parameter:
929 . linesearch - the line search context
930 
931   Output Parameter:
932 . type - The type of line search, or `NULL` if not set
933 
934   Level: intermediate
935 
936 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchSetType()`
937 @*/
938 PetscErrorCode SNESLineSearchGetType(SNESLineSearch linesearch, SNESLineSearchType *type)
939 {
940   PetscFunctionBegin;
941   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
942   PetscAssertPointer(type, 2);
943   *type = ((PetscObject)linesearch)->type_name;
944   PetscFunctionReturn(PETSC_SUCCESS);
945 }
946 
947 /*@C
948   SNESLineSearchSetType - Sets the `SNESLinesearchType` of a `SNESLineSearch`
949 
950   Logically Collective
951 
952   Input Parameters:
953 + linesearch - the line search context
954 - type       - The type of line search to be used, see `SNESLineSearchType`
955 
956   Options Database Key:
957 . -snes_linesearch_type <type> - basic (or equivalently none), bt, l2, cp, nleqerr, shell
958 
959   Level: intermediate
960 
961 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchGetType()`
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 /*@C
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 /*@C
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   Fortran Notes:
1591   The user should pass in a string 'prefix' of
1592   sufficient length to hold the prefix.
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   Input Parameters:
1656 + linesearch - the line search context
1657 - result     - The success or failure status
1658 
1659   Level: developer
1660 
1661   Note:
1662   This is typically called in a `SNESLineSearchType` implementation of `SNESLineSearchApply()` or a `SNESLINESEARCHSHELL` implementation to set
1663   the success or failure of the line search method.
1664 
1665 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchReason`, `SNESLineSearchGetSResult()`
1666 @*/
1667 PetscErrorCode SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result)
1668 {
1669   PetscFunctionBegin;
1670   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1671   linesearch->result = result;
1672   PetscFunctionReturn(PETSC_SUCCESS);
1673 }
1674 
1675 // PetscClangLinter pragma disable: -fdoc-param-list-func-parameter-documentation
1676 /*@C
1677   SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.
1678 
1679   Logically Collective
1680 
1681   Input Parameters:
1682 + linesearch  - the linesearch object
1683 . projectfunc - function for projecting the function to the bounds, see `SNESLineSearchVIProjectFn` for calling sequence
1684 - normfunc    - function for computing the norm of an active set, see `SNESLineSearchVINormFn` for calling sequence
1685 
1686   Level: advanced
1687 
1688   Notes:
1689   The VI solvers require projection of the solution to the feasible set.  `projectfunc` should implement this.
1690 
1691   The VI solvers require special evaluation of the function norm such that the norm is only calculated
1692   on the inactive set.  This should be implemented by `normfunc`.
1693 
1694 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetVIFunctions()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()`,
1695           `SNESLineSearchVIProjectFn`, `SNESLineSearchVINormFn`
1696 @*/
1697 PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFn *projectfunc, SNESLineSearchVINormFn *normfunc)
1698 {
1699   PetscFunctionBegin;
1700   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1701   if (projectfunc) linesearch->ops->viproject = projectfunc;
1702   if (normfunc) linesearch->ops->vinorm = normfunc;
1703   PetscFunctionReturn(PETSC_SUCCESS);
1704 }
1705 
1706 /*@C
1707   SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.
1708 
1709   Not Collective
1710 
1711   Input Parameter:
1712 . linesearch - the line search context, obtain with `SNESGetLineSearch()`
1713 
1714   Output Parameters:
1715 + projectfunc - function for projecting the function to the bounds, see `SNESLineSearchVIProjectFn` for calling sequence
1716 - normfunc    - function for computing the norm of an active set, see `SNESLineSearchVINormFn ` for calling sequence
1717 
1718   Level: advanced
1719 
1720 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetVIFunctions()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()`,
1721           `SNESLineSearchVIProjectFn`, `SNESLineSearchVINormFn`
1722 @*/
1723 PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFn **projectfunc, SNESLineSearchVINormFn **normfunc)
1724 {
1725   PetscFunctionBegin;
1726   if (projectfunc) *projectfunc = linesearch->ops->viproject;
1727   if (normfunc) *normfunc = linesearch->ops->vinorm;
1728   PetscFunctionReturn(PETSC_SUCCESS);
1729 }
1730 
1731 /*@C
1732   SNESLineSearchRegister - register a line search type `SNESLineSearchType`
1733 
1734   Input Parameters:
1735 + sname    - name of the `SNESLineSearchType()`
1736 - function - the creation function for that type
1737 
1738   Calling sequence of `function`:
1739 . ls - the line search context
1740 
1741   Level: advanced
1742 
1743 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchSetType()`
1744 @*/
1745 PetscErrorCode SNESLineSearchRegister(const char sname[], PetscErrorCode (*function)(SNESLineSearch ls))
1746 {
1747   PetscFunctionBegin;
1748   PetscCall(SNESInitializePackage());
1749   PetscCall(PetscFunctionListAdd(&SNESLineSearchList, sname, function));
1750   PetscFunctionReturn(PETSC_SUCCESS);
1751 }
1752