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