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