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