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