xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision 49abdd8a111d9c2ef7fc48ade253ef64e07f9b37)
1af0996ceSBarry Smith #include <petsc/private/linesearchimpl.h> /*I "petscsnes.h" I*/
2bf7f4e0aSPeter Brune 
3f1c6b773SPeter Brune PetscBool         SNESLineSearchRegisterAllCalled = PETSC_FALSE;
40298fd71SBarry Smith PetscFunctionList SNESLineSearchList              = NULL;
5bf7f4e0aSPeter Brune 
6f1c6b773SPeter Brune PetscClassId  SNESLINESEARCH_CLASSID;
757a83faaSBarry Smith PetscLogEvent SNESLINESEARCH_Apply;
8bf7f4e0aSPeter Brune 
9dcf2fd19SBarry Smith /*@
10f6dfbefdSBarry Smith   SNESLineSearchMonitorCancel - Clears all the monitor functions for a `SNESLineSearch` object.
11dcf2fd19SBarry Smith 
12c3339decSBarry Smith   Logically Collective
13dcf2fd19SBarry Smith 
142fe279fdSBarry Smith   Input Parameter:
15f6dfbefdSBarry Smith . ls - the `SNESLineSearch` context
16dcf2fd19SBarry Smith 
17dcf2fd19SBarry Smith   Options Database Key:
18dcf2fd19SBarry Smith . -snes_linesearch_monitor_cancel - cancels all monitors that have been hardwired
19f6dfbefdSBarry Smith     into a code by calls to `SNESLineSearchMonitorSet()`, but does not cancel those
20dcf2fd19SBarry Smith     set via the options database
21dcf2fd19SBarry Smith 
2220f4b53cSBarry Smith   Level: advanced
2320f4b53cSBarry Smith 
24dcf2fd19SBarry Smith   Notes:
25f6dfbefdSBarry Smith   There is no way to clear one specific monitor from a `SNESLineSearch` object.
26dcf2fd19SBarry Smith 
27420bcc1bSBarry Smith   This does not clear the monitor set with `SNESLineSearchSetDefaultMonitor()` use `SNESLineSearchSetDefaultMonitor`(`ls`,`NULL`) to cancel it
28dcf2fd19SBarry Smith   that one.
29dcf2fd19SBarry Smith 
30420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorDefault()`, `SNESLineSearchMonitorSet()`
31dcf2fd19SBarry Smith @*/
32d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchMonitorCancel(SNESLineSearch ls)
33d71ae5a4SJacob Faibussowitsch {
34dcf2fd19SBarry Smith   PetscInt i;
35dcf2fd19SBarry Smith 
36dcf2fd19SBarry Smith   PetscFunctionBegin;
37dcf2fd19SBarry Smith   PetscValidHeaderSpecific(ls, SNESLINESEARCH_CLASSID, 1);
38dcf2fd19SBarry Smith   for (i = 0; i < ls->numbermonitors; i++) {
3948a46eb9SPierre Jolivet     if (ls->monitordestroy[i]) PetscCall((*ls->monitordestroy[i])(&ls->monitorcontext[i]));
40dcf2fd19SBarry Smith   }
41dcf2fd19SBarry Smith   ls->numbermonitors = 0;
423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
43dcf2fd19SBarry Smith }
44dcf2fd19SBarry Smith 
45dcf2fd19SBarry Smith /*@
46dcf2fd19SBarry Smith   SNESLineSearchMonitor - runs the user provided monitor routines, if they exist
47dcf2fd19SBarry Smith 
48c3339decSBarry Smith   Collective
49dcf2fd19SBarry Smith 
502fe279fdSBarry Smith   Input Parameter:
51dcf2fd19SBarry Smith . ls - the linesearch object
52dcf2fd19SBarry Smith 
532fe279fdSBarry Smith   Level: developer
542fe279fdSBarry Smith 
55f6dfbefdSBarry Smith   Note:
56420bcc1bSBarry Smith   This routine is called by the `SNESLineSearch` implementations.
57dcf2fd19SBarry Smith   It does not typically need to be called by the user.
58dcf2fd19SBarry Smith 
59420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorSet()`
60dcf2fd19SBarry Smith @*/
61d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchMonitor(SNESLineSearch ls)
62d71ae5a4SJacob Faibussowitsch {
63dcf2fd19SBarry Smith   PetscInt i, n = ls->numbermonitors;
64dcf2fd19SBarry Smith 
65dcf2fd19SBarry Smith   PetscFunctionBegin;
6648a46eb9SPierre Jolivet   for (i = 0; i < n; i++) PetscCall((*ls->monitorftns[i])(ls, ls->monitorcontext[i]));
673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
68dcf2fd19SBarry Smith }
69dcf2fd19SBarry Smith 
70dcf2fd19SBarry Smith /*@C
71dcf2fd19SBarry Smith   SNESLineSearchMonitorSet - Sets an ADDITIONAL function that is to be used at every
72dcf2fd19SBarry Smith   iteration of the nonlinear solver to display the iteration's
73dcf2fd19SBarry Smith   progress.
74dcf2fd19SBarry Smith 
75c3339decSBarry Smith   Logically Collective
76dcf2fd19SBarry Smith 
77dcf2fd19SBarry Smith   Input Parameters:
78f6dfbefdSBarry Smith + ls             - the `SNESLineSearch` context
79dcf2fd19SBarry Smith . f              - the monitor function
80420bcc1bSBarry Smith . mctx           - [optional] user-defined context for private data for the monitor routine (use `NULL` if no context is desired)
81*49abdd8aSBarry Smith - monitordestroy - [optional] routine that frees monitor context (may be `NULL`), see `PetscCtxDestroyFn` for the calling sequence
82420bcc1bSBarry Smith 
83420bcc1bSBarry Smith   Calling sequence of `f`:
84420bcc1bSBarry Smith + ls   - the `SNESLineSearch` context
85420bcc1bSBarry Smith - mctx - [optional] user-defined context for private data for the monitor routine
86420bcc1bSBarry Smith 
8720f4b53cSBarry Smith   Level: intermediate
88dcf2fd19SBarry Smith 
89f6dfbefdSBarry Smith   Note:
90dcf2fd19SBarry Smith   Several different monitoring routines may be set by calling
91f6dfbefdSBarry Smith   `SNESLineSearchMonitorSet()` multiple times; all will be called in the
92dcf2fd19SBarry Smith   order in which they were set.
93dcf2fd19SBarry Smith 
94420bcc1bSBarry Smith   Fortran Note:
95f6dfbefdSBarry Smith   Only a single monitor function can be set for each `SNESLineSearch` object
96dcf2fd19SBarry Smith 
97*49abdd8aSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorDefault()`, `SNESLineSearchMonitorCancel()`, `PetscCtxDestroyFn`
98dcf2fd19SBarry Smith @*/
99*49abdd8aSBarry Smith PetscErrorCode SNESLineSearchMonitorSet(SNESLineSearch ls, PetscErrorCode (*f)(SNESLineSearch ls, void *mctx), void *mctx, PetscCtxDestroyFn *monitordestroy)
100d71ae5a4SJacob Faibussowitsch {
10178064530SBarry Smith   PetscInt  i;
10278064530SBarry Smith   PetscBool identical;
10378064530SBarry Smith 
104dcf2fd19SBarry Smith   PetscFunctionBegin;
105dcf2fd19SBarry Smith   PetscValidHeaderSpecific(ls, SNESLINESEARCH_CLASSID, 1);
10678064530SBarry Smith   for (i = 0; i < ls->numbermonitors; i++) {
1079566063dSJacob Faibussowitsch     PetscCall(PetscMonitorCompare((PetscErrorCode (*)(void))f, mctx, monitordestroy, (PetscErrorCode (*)(void))ls->monitorftns[i], ls->monitorcontext[i], ls->monitordestroy[i], &identical));
1083ba16761SJacob Faibussowitsch     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
10978064530SBarry Smith   }
1105f80ce2aSJacob Faibussowitsch   PetscCheck(ls->numbermonitors < MAXSNESLSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many monitors set");
111dcf2fd19SBarry Smith   ls->monitorftns[ls->numbermonitors]      = f;
112dcf2fd19SBarry Smith   ls->monitordestroy[ls->numbermonitors]   = monitordestroy;
113dcf2fd19SBarry Smith   ls->monitorcontext[ls->numbermonitors++] = (void *)mctx;
1143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
115dcf2fd19SBarry Smith }
116dcf2fd19SBarry Smith 
117dcf2fd19SBarry Smith /*@C
118f6dfbefdSBarry Smith   SNESLineSearchMonitorSolutionUpdate - Monitors each update of the function value the linesearch tries
119dcf2fd19SBarry Smith 
120c3339decSBarry Smith   Collective
121dcf2fd19SBarry Smith 
122dcf2fd19SBarry Smith   Input Parameters:
123420bcc1bSBarry Smith + ls - the `SNESLineSearch` object
124f6dfbefdSBarry Smith - vf - the context for the monitor, in this case it is an `PetscViewerAndFormat`
125dcf2fd19SBarry Smith 
126420bcc1bSBarry Smith   Options Database Key:
127420bcc1bSBarry Smith . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine
128dcf2fd19SBarry Smith 
129420bcc1bSBarry Smith   Level: developer
130420bcc1bSBarry Smith 
131420bcc1bSBarry Smith   This is not normally called directly but is passed to `SNESLineSearchMonitorSet()`
132420bcc1bSBarry Smith 
133420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorSet()`, `SNESMonitorSolution()`
134dcf2fd19SBarry Smith @*/
135d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchMonitorSolutionUpdate(SNESLineSearch ls, PetscViewerAndFormat *vf)
136d71ae5a4SJacob Faibussowitsch {
137d12e167eSBarry Smith   PetscViewer viewer = vf->viewer;
138dcf2fd19SBarry Smith   Vec         Y, W, G;
139dcf2fd19SBarry Smith 
140dcf2fd19SBarry Smith   PetscFunctionBegin;
1419566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetVecs(ls, NULL, NULL, &Y, &W, &G));
1429566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(viewer, vf->format));
1439566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted update to solution \n"));
1449566063dSJacob Faibussowitsch   PetscCall(VecView(Y, viewer));
1459566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted new solution \n"));
1469566063dSJacob Faibussowitsch   PetscCall(VecView(W, viewer));
1479566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted updated function value\n"));
1489566063dSJacob Faibussowitsch   PetscCall(VecView(G, viewer));
1499566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(viewer));
1503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
151dcf2fd19SBarry Smith }
152dcf2fd19SBarry Smith 
153f40b411bSPeter Brune /*@
154420bcc1bSBarry Smith   SNESLineSearchCreate - Creates a `SNESLineSearch` context.
155f40b411bSPeter Brune 
156f6dfbefdSBarry Smith   Logically Collective
157f40b411bSPeter Brune 
1582fe279fdSBarry Smith   Input Parameter:
159f6dfbefdSBarry Smith . comm - MPI communicator for the line search (typically from the associated `SNES` context).
160f40b411bSPeter Brune 
1612fe279fdSBarry Smith   Output Parameter:
1628e557f58SPeter Brune . outlinesearch - the new line search context
163f40b411bSPeter Brune 
164162e0bf5SPeter Brune   Level: developer
165162e0bf5SPeter Brune 
166f6dfbefdSBarry Smith   Note:
167420bcc1bSBarry Smith   The preferred calling sequence is to use `SNESGetLineSearch()` to acquire the `SNESLineSearch` instance
168f6dfbefdSBarry Smith   already associated with the `SNES`.
169f40b411bSPeter Brune 
170420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `LineSearchDestroy()`, `SNESGetLineSearch()`
171f40b411bSPeter Brune @*/
172d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch)
173d71ae5a4SJacob Faibussowitsch {
174f1c6b773SPeter Brune   SNESLineSearch linesearch;
175bf388a1fSBarry Smith 
176bf7f4e0aSPeter Brune   PetscFunctionBegin;
1774f572ea9SToby Isaac   PetscAssertPointer(outlinesearch, 2);
1789566063dSJacob Faibussowitsch   PetscCall(SNESInitializePackage());
179f5af7f23SKarl Rupp 
1809566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(linesearch, SNESLINESEARCH_CLASSID, "SNESLineSearch", "Linesearch", "SNESLineSearch", comm, SNESLineSearchDestroy, SNESLineSearchView));
1810298fd71SBarry Smith   linesearch->vec_sol_new  = NULL;
1820298fd71SBarry Smith   linesearch->vec_func_new = NULL;
1830298fd71SBarry Smith   linesearch->vec_sol      = NULL;
1840298fd71SBarry Smith   linesearch->vec_func     = NULL;
1850298fd71SBarry Smith   linesearch->vec_update   = NULL;
1869bd66eb0SPeter Brune 
187bf7f4e0aSPeter Brune   linesearch->lambda       = 1.0;
188bf7f4e0aSPeter Brune   linesearch->fnorm        = 1.0;
189bf7f4e0aSPeter Brune   linesearch->ynorm        = 1.0;
190bf7f4e0aSPeter Brune   linesearch->xnorm        = 1.0;
191422a814eSBarry Smith   linesearch->result       = SNES_LINESEARCH_SUCCEEDED;
192bf7f4e0aSPeter Brune   linesearch->norms        = PETSC_TRUE;
193bf7f4e0aSPeter Brune   linesearch->keeplambda   = PETSC_FALSE;
194bf7f4e0aSPeter Brune   linesearch->damping      = 1.0;
195bf7f4e0aSPeter Brune   linesearch->maxstep      = 1e8;
196bf7f4e0aSPeter Brune   linesearch->steptol      = 1e-12;
197516fe3c3SPeter Brune   linesearch->rtol         = 1e-8;
198516fe3c3SPeter Brune   linesearch->atol         = 1e-15;
199516fe3c3SPeter Brune   linesearch->ltol         = 1e-8;
2000298fd71SBarry Smith   linesearch->precheckctx  = NULL;
2010298fd71SBarry Smith   linesearch->postcheckctx = NULL;
20259405d5eSPeter Brune   linesearch->max_its      = 1;
203bf7f4e0aSPeter Brune   linesearch->setupcalled  = PETSC_FALSE;
2043add74b1SBarry Smith   linesearch->monitor      = NULL;
205bf7f4e0aSPeter Brune   *outlinesearch           = linesearch;
2063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
207bf7f4e0aSPeter Brune }
208bf7f4e0aSPeter Brune 
209f40b411bSPeter Brune /*@
21078bcb3b5SPeter Brune   SNESLineSearchSetUp - Prepares the line search for being applied by allocating
21178bcb3b5SPeter Brune   any required vectors.
212f40b411bSPeter Brune 
213c3339decSBarry Smith   Collective
214f40b411bSPeter Brune 
2152fe279fdSBarry Smith   Input Parameter:
216f6dfbefdSBarry Smith . linesearch - The `SNESLineSearch` instance.
217f40b411bSPeter Brune 
21820f4b53cSBarry Smith   Level: advanced
21920f4b53cSBarry Smith 
220f6dfbefdSBarry Smith   Note:
221f6dfbefdSBarry Smith   For most cases, this needn't be called by users or outside of `SNESLineSearchApply()`.
222cd7522eaSPeter Brune   The only current case where this is called outside of this is for the VI
22378bcb3b5SPeter Brune   solvers, which modify the solution and work vectors before the first call
224f6dfbefdSBarry Smith   of `SNESLineSearchApply()`, requiring the `SNESLineSearch` work vectors to be
225cd7522eaSPeter Brune   allocated upfront.
226cd7522eaSPeter Brune 
227420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchReset()`
228f40b411bSPeter Brune @*/
229d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch)
230d71ae5a4SJacob Faibussowitsch {
231bf7f4e0aSPeter Brune   PetscFunctionBegin;
23248a46eb9SPierre Jolivet   if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC));
233bf7f4e0aSPeter Brune   if (!linesearch->setupcalled) {
23448a46eb9SPierre Jolivet     if (!linesearch->vec_sol_new) PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new));
23548a46eb9SPierre Jolivet     if (!linesearch->vec_func_new) PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new));
236213acdd3SPierre Jolivet     PetscTryTypeMethod(linesearch, setup);
2379566063dSJacob Faibussowitsch     if (!linesearch->ops->snesfunc) PetscCall(SNESLineSearchSetFunction(linesearch, SNESComputeFunction));
238bf7f4e0aSPeter Brune     linesearch->lambda      = linesearch->damping;
239bf7f4e0aSPeter Brune     linesearch->setupcalled = PETSC_TRUE;
240bf7f4e0aSPeter Brune   }
2413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
242bf7f4e0aSPeter Brune }
243bf7f4e0aSPeter Brune 
244f40b411bSPeter Brune /*@
245f6dfbefdSBarry Smith   SNESLineSearchReset - Undoes the `SNESLineSearchSetUp()` and deletes any `Vec`s or `Mat`s allocated by the line search.
246f40b411bSPeter Brune 
247c3339decSBarry Smith   Collective
248f40b411bSPeter Brune 
2492fe279fdSBarry Smith   Input Parameter:
250f6dfbefdSBarry Smith . linesearch - The `SNESLineSearch` instance.
251f40b411bSPeter Brune 
25220f4b53cSBarry Smith   Level: developer
25320f4b53cSBarry Smith 
254f6dfbefdSBarry Smith   Note:
255f6dfbefdSBarry Smith   Usually only called by `SNESReset()`
256f190f2fcSBarry Smith 
257420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchSetUp()`
258f40b411bSPeter Brune @*/
259d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch)
260d71ae5a4SJacob Faibussowitsch {
261bf7f4e0aSPeter Brune   PetscFunctionBegin;
262213acdd3SPierre Jolivet   PetscTryTypeMethod(linesearch, reset);
263f5af7f23SKarl Rupp 
2649566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&linesearch->vec_sol_new));
2659566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&linesearch->vec_func_new));
266bf7f4e0aSPeter Brune 
2679566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(linesearch->nwork, &linesearch->work));
268f5af7f23SKarl Rupp 
269bf7f4e0aSPeter Brune   linesearch->nwork       = 0;
270bf7f4e0aSPeter Brune   linesearch->setupcalled = PETSC_FALSE;
2713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
272bf7f4e0aSPeter Brune }
273bf7f4e0aSPeter Brune 
274ed07d7d7SPeter Brune /*@C
275f6dfbefdSBarry Smith   SNESLineSearchSetFunction - Sets the function evaluation used by the `SNES` line search
276f6dfbefdSBarry Smith   `
277e4094ef1SJacob Faibussowitsch 
278ed07d7d7SPeter Brune   Input Parameters:
279e4094ef1SJacob Faibussowitsch + linesearch - the `SNESLineSearch` context
280e4094ef1SJacob Faibussowitsch - func       - function evaluation routine, this is usually the function provided with `SNESSetFunction()`
281ed07d7d7SPeter Brune 
282420bcc1bSBarry Smith   Calling sequence of `func`:
283420bcc1bSBarry Smith + snes - the `SNES` with which the `SNESLineSearch` context is associated with
284420bcc1bSBarry Smith . x    - the input vector
285420bcc1bSBarry Smith - f    - the computed value of the function
286420bcc1bSBarry Smith 
287ed07d7d7SPeter Brune   Level: developer
288ed07d7d7SPeter Brune 
289420bcc1bSBarry Smith   Note:
290420bcc1bSBarry Smith   By default the `SNESLineSearch` uses the function provided by `SNESSetFunction()` so this is rarely needed
291420bcc1bSBarry Smith 
292420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESSetFunction()`
293ed07d7d7SPeter Brune @*/
294420bcc1bSBarry Smith PetscErrorCode SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES snes, Vec x, Vec f))
295d71ae5a4SJacob Faibussowitsch {
296ed07d7d7SPeter Brune   PetscFunctionBegin;
297ed07d7d7SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
298ed07d7d7SPeter Brune   linesearch->ops->snesfunc = func;
2993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
300ed07d7d7SPeter Brune }
301ed07d7d7SPeter Brune 
30286d74e61SPeter Brune /*@C
303420bcc1bSBarry Smith   SNESLineSearchSetPreCheck - Sets a function that is called after the initial search direction has been computed but
304420bcc1bSBarry Smith   before the line search routine has been applied. Allows adjusting the result of (usually a linear solve) that
305f190f2fcSBarry Smith   determined the search direction.
30686d74e61SPeter Brune 
307c3339decSBarry Smith   Logically Collective
30886d74e61SPeter Brune 
30986d74e61SPeter Brune   Input Parameters:
310f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context
311420bcc1bSBarry Smith . func       - [optional] function evaluation routine
31220f4b53cSBarry Smith - ctx        - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
31386d74e61SPeter Brune 
314420bcc1bSBarry Smith   Calling sequence of `func`:
315420bcc1bSBarry Smith + ls        - the `SNESLineSearch` context
316420bcc1bSBarry Smith . x         - the current solution
317a678f235SPierre Jolivet . d         - the current search direction
318420bcc1bSBarry Smith . changed_d - indicates if the search direction has been changed
319420bcc1bSBarry Smith - ctx       - the context passed to `SNESLineSearchSetPreCheck()`
320420bcc1bSBarry Smith 
32186d74e61SPeter Brune   Level: intermediate
32286d74e61SPeter Brune 
323f6dfbefdSBarry Smith   Note:
324420bcc1bSBarry Smith   Use `SNESLineSearchSetPostCheck()` to change the step after the line search is complete.
325f0b84518SBarry Smith 
326f0b84518SBarry Smith   Use `SNESVISetVariableBounds()` and `SNESVISetComputeVariableBounds()` to cause `SNES` to automatically control the ranges of variables allowed.
327f0b84518SBarry Smith 
328420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()`,
329869ba2dcSBarry Smith           `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`
330f0b84518SBarry Smith 
33186d74e61SPeter Brune @*/
332420bcc1bSBarry Smith PetscErrorCode SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch ls, Vec x, Vec d, PetscBool *changed_d, void *ctx), void *ctx)
333d71ae5a4SJacob Faibussowitsch {
3349bd66eb0SPeter Brune   PetscFunctionBegin;
335f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
336f190f2fcSBarry Smith   if (func) linesearch->ops->precheck = func;
33786d74e61SPeter Brune   if (ctx) linesearch->precheckctx = ctx;
3383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
33986d74e61SPeter Brune }
34086d74e61SPeter Brune 
34186d74e61SPeter Brune /*@C
34253deb899SBarry Smith   SNESLineSearchGetPreCheck - Gets the pre-check function for the line search routine.
34386d74e61SPeter Brune 
344f899ff85SJose E. Roman   Input Parameter:
345f6dfbefdSBarry Smith . linesearch - the `SNESLineSearch` context
34686d74e61SPeter Brune 
34786d74e61SPeter Brune   Output Parameters:
348420bcc1bSBarry Smith + func - [optional] function evaluation routine,  for calling sequence see `SNESLineSearchSetPreCheck()`
34920f4b53cSBarry Smith - ctx  - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
35086d74e61SPeter Brune 
35186d74e61SPeter Brune   Level: intermediate
35286d74e61SPeter Brune 
353420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchSetPostCheck()`
35486d74e61SPeter Brune @*/
355d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch, Vec, Vec, PetscBool *, void *), void **ctx)
356d71ae5a4SJacob Faibussowitsch {
3579bd66eb0SPeter Brune   PetscFunctionBegin;
358f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
359f190f2fcSBarry Smith   if (func) *func = linesearch->ops->precheck;
36086d74e61SPeter Brune   if (ctx) *ctx = linesearch->precheckctx;
3613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
36286d74e61SPeter Brune }
36386d74e61SPeter Brune 
36486d74e61SPeter Brune /*@C
365f190f2fcSBarry Smith   SNESLineSearchSetPostCheck - Sets a user function that is called after the line search has been applied to determine the step
366f190f2fcSBarry Smith   direction and length. Allows the user a chance to change or override the decision of the line search routine
36786d74e61SPeter Brune 
36820f4b53cSBarry Smith   Logically Collective
36986d74e61SPeter Brune 
37086d74e61SPeter Brune   Input Parameters:
371f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context
372420bcc1bSBarry Smith . func       - [optional] function evaluation routine
37320f4b53cSBarry Smith - ctx        - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
37486d74e61SPeter Brune 
375420bcc1bSBarry Smith   Calling sequence of `func`:
376420bcc1bSBarry Smith + ls        - the `SNESLineSearch` context
377420bcc1bSBarry Smith . x         - the current solution
378a678f235SPierre Jolivet . d         - the current search direction
379a678f235SPierre Jolivet . w         - $ w = x + lambda*d $ for some lambda
380420bcc1bSBarry Smith . changed_d - indicates if the search direction `d` has been changed
381420bcc1bSBarry Smith . changed_w - indicates `w` has been changed
382420bcc1bSBarry Smith - ctx       - the context passed to `SNESLineSearchSetPreCheck()`
383420bcc1bSBarry Smith 
38486d74e61SPeter Brune   Level: intermediate
38586d74e61SPeter Brune 
386f0b84518SBarry Smith   Notes:
3876b095a85SStefano Zampini   Use `SNESLineSearchSetPreCheck()` to change the step before the line search is completed.
3886b095a85SStefano Zampini   The calling sequence of the callback does not contain the current scaling factor. To access the value, use `SNESLineSearchGetLambda()`.
389f0b84518SBarry Smith 
390f0b84518SBarry Smith   Use `SNESVISetVariableBounds()` and `SNESVISetComputeVariableBounds()` to cause `SNES` to automatically control the ranges of variables allowed.
391f0b84518SBarry Smith 
392420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchGetPostCheck()`,
393e4094ef1SJacob Faibussowitsch           `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`
39486d74e61SPeter Brune @*/
395420bcc1bSBarry Smith PetscErrorCode SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch ls, Vec x, Vec d, Vec w, PetscBool *changed_d, PetscBool *changed_w, void *ctx), void *ctx)
396d71ae5a4SJacob Faibussowitsch {
39786d74e61SPeter Brune   PetscFunctionBegin;
398f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
399f190f2fcSBarry Smith   if (func) linesearch->ops->postcheck = func;
40086d74e61SPeter Brune   if (ctx) linesearch->postcheckctx = ctx;
4013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
40286d74e61SPeter Brune }
40386d74e61SPeter Brune 
40486d74e61SPeter Brune /*@C
405f1c6b773SPeter Brune   SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine.
40686d74e61SPeter Brune 
407f899ff85SJose E. Roman   Input Parameter:
408f6dfbefdSBarry Smith . linesearch - the `SNESLineSearch` context
40986d74e61SPeter Brune 
41086d74e61SPeter Brune   Output Parameters:
411420bcc1bSBarry Smith + func - [optional] function evaluation routine, see for the calling sequence `SNESLineSearchSetPostCheck()`
41220f4b53cSBarry Smith - ctx  - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
41386d74e61SPeter Brune 
41486d74e61SPeter Brune   Level: intermediate
41586d74e61SPeter Brune 
416420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`
41786d74e61SPeter Brune @*/
418d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx)
419d71ae5a4SJacob Faibussowitsch {
4209bd66eb0SPeter Brune   PetscFunctionBegin;
421f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
422f190f2fcSBarry Smith   if (func) *func = linesearch->ops->postcheck;
42386d74e61SPeter Brune   if (ctx) *ctx = linesearch->postcheckctx;
4243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
42586d74e61SPeter Brune }
42686d74e61SPeter Brune 
427f40b411bSPeter Brune /*@
428f1c6b773SPeter Brune   SNESLineSearchPreCheck - Prepares the line search for being applied.
429f40b411bSPeter Brune 
430c3339decSBarry Smith   Logically Collective
431f40b411bSPeter Brune 
432f40b411bSPeter Brune   Input Parameters:
4337b1df9c1SPeter Brune + linesearch - The linesearch instance.
4347b1df9c1SPeter Brune . X          - The current solution
4357b1df9c1SPeter Brune - Y          - The step direction
436f40b411bSPeter Brune 
4372fe279fdSBarry Smith   Output Parameter:
438420bcc1bSBarry Smith . changed - Indicator that the precheck routine has changed `Y`
439f40b411bSPeter Brune 
440f0b84518SBarry Smith   Level: advanced
441f40b411bSPeter Brune 
442f6dfbefdSBarry Smith   Note:
443420bcc1bSBarry Smith   This calls any function provided with `SNESLineSearchSetPreCheck()` and is called automatically inside the line search routines
444f6dfbefdSBarry Smith 
445420bcc1bSBarry Smith   Developer Note:
446420bcc1bSBarry Smith   The use of `PetscObjectGetState()` would eliminate the need for the `changed` argument to be provided
447420bcc1bSBarry Smith 
448420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`,
449fe8e7dddSPierre Jolivet           `SNESLineSearchGetPostCheck()`
450f40b411bSPeter Brune @*/
451d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed)
452d71ae5a4SJacob Faibussowitsch {
453bf7f4e0aSPeter Brune   PetscFunctionBegin;
454bf7f4e0aSPeter Brune   *changed = PETSC_FALSE;
4556b2b7091SBarry Smith   if (linesearch->ops->precheck) {
456dbbe0bcdSBarry Smith     PetscUseTypeMethod(linesearch, precheck, X, Y, changed, linesearch->precheckctx);
45738bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch, *changed, 4);
458bf7f4e0aSPeter Brune   }
4593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
460bf7f4e0aSPeter Brune }
461bf7f4e0aSPeter Brune 
462f40b411bSPeter Brune /*@
463ef46b1a6SStefano Zampini   SNESLineSearchPostCheck - Hook to modify step direction or updated solution after a successful linesearch
464f40b411bSPeter Brune 
465c3339decSBarry Smith   Logically Collective
466f40b411bSPeter Brune 
467f40b411bSPeter Brune   Input Parameters:
4687b1df9c1SPeter Brune + linesearch - The line search context
4697b1df9c1SPeter Brune . X          - The last solution
4707b1df9c1SPeter Brune . Y          - The step direction
4716b095a85SStefano Zampini - W          - The updated solution, `W = X - lambda * Y` for some lambda
472f40b411bSPeter Brune 
473f40b411bSPeter Brune   Output Parameters:
4746b095a85SStefano Zampini + changed_Y - Indicator if the direction `Y` has been changed.
475420bcc1bSBarry Smith - changed_W - Indicator if the new candidate solution `W` has been changed.
476f40b411bSPeter Brune 
47720f4b53cSBarry Smith   Level: developer
47820f4b53cSBarry Smith 
479f6dfbefdSBarry Smith   Note:
480420bcc1bSBarry Smith   This calls any function provided with `SNESLineSearchSetPostCheck()` and is called automatically inside the line search routines
481f6dfbefdSBarry Smith 
482420bcc1bSBarry Smith   Developer Note:
483420bcc1bSBarry Smith   The use of `PetscObjectGetState()` would eliminate the need for the `changed_Y` and `changed_W` arguments to be provided
484420bcc1bSBarry Smith 
485420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPrecheck()`, `SNESLineSearchGetPrecheck()`
486f40b411bSPeter Brune @*/
487d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W)
488d71ae5a4SJacob Faibussowitsch {
489bf7f4e0aSPeter Brune   PetscFunctionBegin;
490bf7f4e0aSPeter Brune   *changed_Y = PETSC_FALSE;
491bf7f4e0aSPeter Brune   *changed_W = PETSC_FALSE;
4926b2b7091SBarry Smith   if (linesearch->ops->postcheck) {
493dbbe0bcdSBarry Smith     PetscUseTypeMethod(linesearch, postcheck, X, Y, W, changed_Y, changed_W, linesearch->postcheckctx);
49438bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch, *changed_Y, 5);
49538bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch, *changed_W, 6);
49686d74e61SPeter Brune   }
4973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
49886d74e61SPeter Brune }
49986d74e61SPeter Brune 
50086d74e61SPeter Brune /*@C
5011d27aa22SBarry Smith   SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration {cite}`hindmarsh1996time`
50286d74e61SPeter Brune 
503c3339decSBarry Smith   Logically Collective
50486d74e61SPeter Brune 
5054165533cSJose E. Roman   Input Parameters:
506420bcc1bSBarry Smith + linesearch - the line search context
50786d74e61SPeter Brune . X          - base state for this step
508907376e6SBarry Smith - ctx        - context for this function
50986d74e61SPeter Brune 
51097bb3fdcSJose E. Roman   Input/Output Parameter:
51197bb3fdcSJose E. Roman . Y - correction, possibly modified
51297bb3fdcSJose E. Roman 
51397bb3fdcSJose E. Roman   Output Parameter:
514420bcc1bSBarry Smith . changed - flag indicating that `Y` was modified
51586d74e61SPeter Brune 
516420bcc1bSBarry Smith   Options Database Keys:
517cd7522eaSPeter Brune + -snes_linesearch_precheck_picard       - activate this routine
518cd7522eaSPeter Brune - -snes_linesearch_precheck_picard_angle - angle
51986d74e61SPeter Brune 
52086d74e61SPeter Brune   Level: advanced
52186d74e61SPeter Brune 
52286d74e61SPeter Brune   Notes:
523f6dfbefdSBarry Smith   This function should be passed to `SNESLineSearchSetPreCheck()`
52486d74e61SPeter Brune 
52586d74e61SPeter Brune   The justification for this method involves the linear convergence of a Picard iteration
5261d27aa22SBarry Smith   so the Picard linearization should be provided in place of the "Jacobian"  {cite}`hindmarsh1996time`. This correction
52786d74e61SPeter Brune   is generally not useful when using a Newton linearization.
52886d74e61SPeter Brune 
529420bcc1bSBarry Smith   Developer Note:
530420bcc1bSBarry Smith   The use of `PetscObjectGetState()` would eliminate the need for the `changed` argument to be provided
531420bcc1bSBarry Smith 
532420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESSetPicard()`, `SNESGetLineSearch()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchSetPostCheck()`
53386d74e61SPeter Brune @*/
534d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed, void *ctx)
535d71ae5a4SJacob Faibussowitsch {
53686d74e61SPeter Brune   PetscReal   angle = *(PetscReal *)linesearch->precheckctx;
53786d74e61SPeter Brune   Vec         Ylast;
53886d74e61SPeter Brune   PetscScalar dot;
53986d74e61SPeter Brune   PetscInt    iter;
54086d74e61SPeter Brune   PetscReal   ynorm, ylastnorm, theta, angle_radians;
54186d74e61SPeter Brune   SNES        snes;
54286d74e61SPeter Brune 
54386d74e61SPeter Brune   PetscFunctionBegin;
5449566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
5459566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)snes, "SNESLineSearchPreCheckPicard_Ylast", (PetscObject *)&Ylast));
54686d74e61SPeter Brune   if (!Ylast) {
5479566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(Y, &Ylast));
5489566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)snes, "SNESLineSearchPreCheckPicard_Ylast", (PetscObject)Ylast));
5499566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)Ylast));
55086d74e61SPeter Brune   }
5519566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &iter));
55286d74e61SPeter Brune   if (iter < 2) {
5539566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y, Ylast));
55486d74e61SPeter Brune     *changed = PETSC_FALSE;
5553ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
55686d74e61SPeter Brune   }
55786d74e61SPeter Brune 
5589566063dSJacob Faibussowitsch   PetscCall(VecDot(Y, Ylast, &dot));
5599566063dSJacob Faibussowitsch   PetscCall(VecNorm(Y, NORM_2, &ynorm));
5609566063dSJacob Faibussowitsch   PetscCall(VecNorm(Ylast, NORM_2, &ylastnorm));
56186d74e61SPeter Brune   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
562255453a1SBarry Smith   theta         = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm), -1.0, 1.0));
56386d74e61SPeter Brune   angle_radians = angle * PETSC_PI / 180.;
56486d74e61SPeter Brune   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
56586d74e61SPeter Brune     /* Modify the step Y */
56686d74e61SPeter Brune     PetscReal alpha, ydiffnorm;
5679566063dSJacob Faibussowitsch     PetscCall(VecAXPY(Ylast, -1.0, Y));
5689566063dSJacob Faibussowitsch     PetscCall(VecNorm(Ylast, NORM_2, &ydiffnorm));
569f85e2ce2SBarry Smith     alpha = (ydiffnorm > .001 * ylastnorm) ? ylastnorm / ydiffnorm : 1000.0;
5709566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y, Ylast));
5719566063dSJacob Faibussowitsch     PetscCall(VecScale(Y, alpha));
5729566063dSJacob Faibussowitsch     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));
573a47ec85fSBarry Smith     *changed = PETSC_TRUE;
57486d74e61SPeter Brune   } else {
5759566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Angle %14.12e degrees exceeds threshold %14.12e, no correction applied\n", (double)(theta * 180. / PETSC_PI), (double)angle));
5769566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y, Ylast));
57786d74e61SPeter Brune     *changed = PETSC_FALSE;
578bf7f4e0aSPeter Brune   }
5793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
580bf7f4e0aSPeter Brune }
581bf7f4e0aSPeter Brune 
582f40b411bSPeter Brune /*@
583cd7522eaSPeter Brune   SNESLineSearchApply - Computes the line-search update.
584f40b411bSPeter Brune 
585c3339decSBarry Smith   Collective
586f40b411bSPeter Brune 
5876b095a85SStefano Zampini   Input Parameter:
5886b095a85SStefano Zampini . linesearch - The line search context
589f40b411bSPeter Brune 
5906b867d5aSJose E. Roman   Input/Output Parameters:
5916b867d5aSJose E. Roman + X     - The current solution, on output the new solution
592420bcc1bSBarry Smith . F     - The current function value, on output the new function value at the solution value `X`
5931bd63e3eSSatish Balay . fnorm - The current norm of `F`, on output the new norm of `F`
5946b095a85SStefano Zampini - Y     - The current search direction, on output the direction determined by the linesearch, i.e. Xnew = Xold - lambda*Y
595f40b411bSPeter Brune 
596cd7522eaSPeter Brune   Options Database Keys:
5970b00b554SBarry Smith + -snes_linesearch_type                - basic (or equivalently none), bt, l2, cp, nleqerr, shell
598dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches
5991fe24845SBarry Smith . -snes_linesearch_damping             - The linesearch damping parameter, default is 1.0 (no damping)
6001fe24845SBarry Smith . -snes_linesearch_norms               - Turn on/off the linesearch norms computation (SNESLineSearchSetComputeNorms())
6013c7d6663SPeter Brune . -snes_linesearch_keeplambda          - Keep the previous search length as the initial guess
6023c7d6663SPeter Brune - -snes_linesearch_max_it              - The number of iterations for iterative line searches
603cd7522eaSPeter Brune 
604e4094ef1SJacob Faibussowitsch   Level: intermediate
60520f4b53cSBarry Smith 
606cd7522eaSPeter Brune   Notes:
607f6dfbefdSBarry Smith   This is typically called from within a `SNESSolve()` implementation in order to
608f6dfbefdSBarry Smith   help with convergence of the nonlinear method.  Various `SNES` types use line searches
609cd7522eaSPeter Brune   in different ways, but the overarching theme is that a line search is used to determine
610cd7522eaSPeter Brune   an optimal damping parameter of a step at each iteration of the method.  Each
611f6dfbefdSBarry Smith   application of the line search may invoke `SNESComputeFunction()` several times, and
612cd7522eaSPeter Brune   therefore may be fairly expensive.
613cd7522eaSPeter Brune 
6141bd63e3eSSatish Balay .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchGetLambda()`, `SNESLineSearchPreCheck()`, `SNESLineSearchPostCheck()`, `SNESSolve()`, `SNESComputeFunction()`, `SNESLineSearchSetComputeNorms()`,
615db781477SPatrick Sanan           `SNESLineSearchType`, `SNESLineSearchSetType()`
616f40b411bSPeter Brune @*/
617d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal *fnorm, Vec Y)
618d71ae5a4SJacob Faibussowitsch {
619bf388a1fSBarry Smith   PetscFunctionBegin;
620f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
621bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
622bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
623064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(Y, VEC_CLASSID, 5);
624bf7f4e0aSPeter Brune 
625422a814eSBarry Smith   linesearch->result = SNES_LINESEARCH_SUCCEEDED;
626bf7f4e0aSPeter Brune 
627bf7f4e0aSPeter Brune   linesearch->vec_sol    = X;
628bf7f4e0aSPeter Brune   linesearch->vec_update = Y;
629bf7f4e0aSPeter Brune   linesearch->vec_func   = F;
630bf7f4e0aSPeter Brune 
6319566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchSetUp(linesearch));
632bf7f4e0aSPeter Brune 
633f5af7f23SKarl Rupp   if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */
634bf7f4e0aSPeter Brune 
635f5af7f23SKarl Rupp   if (fnorm) linesearch->fnorm = *fnorm;
63648a46eb9SPierre Jolivet   else PetscCall(VecNorm(F, NORM_2, &linesearch->fnorm));
637bf7f4e0aSPeter Brune 
6389566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(SNESLINESEARCH_Apply, linesearch, X, F, Y));
639bf7f4e0aSPeter Brune 
640dbbe0bcdSBarry Smith   PetscUseTypeMethod(linesearch, apply);
641bf7f4e0aSPeter Brune 
6429566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(SNESLINESEARCH_Apply, linesearch, X, F, Y));
643bf7f4e0aSPeter Brune 
644f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
6453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
646bf7f4e0aSPeter Brune }
647bf7f4e0aSPeter Brune 
648f40b411bSPeter Brune /*@
649f1c6b773SPeter Brune   SNESLineSearchDestroy - Destroys the line search instance.
650f40b411bSPeter Brune 
651c3339decSBarry Smith   Collective
652f40b411bSPeter Brune 
6532fe279fdSBarry Smith   Input Parameter:
6548e557f58SPeter Brune . linesearch - The line search context
655f40b411bSPeter Brune 
65684238204SBarry Smith   Level: developer
657f40b411bSPeter Brune 
658420bcc1bSBarry Smith   Note:
659420bcc1bSBarry Smith   The line search in `SNES` is automatically called on `SNESDestroy()` so this call is rarely needed
660420bcc1bSBarry Smith 
661420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchReset()`, `SNESDestroy()`
662f40b411bSPeter Brune @*/
663d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchDestroy(SNESLineSearch *linesearch)
664d71ae5a4SJacob Faibussowitsch {
665bf7f4e0aSPeter Brune   PetscFunctionBegin;
6663ba16761SJacob Faibussowitsch   if (!*linesearch) PetscFunctionReturn(PETSC_SUCCESS);
667f4f49eeaSPierre Jolivet   PetscValidHeaderSpecific(*linesearch, SNESLINESEARCH_CLASSID, 1);
668f4f49eeaSPierre Jolivet   if (--((PetscObject)*linesearch)->refct > 0) {
6699371c9d4SSatish Balay     *linesearch = NULL;
6703ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
6719371c9d4SSatish Balay   }
6729566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsViewOff((PetscObject)*linesearch));
6739566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchReset(*linesearch));
6743ba16761SJacob Faibussowitsch   PetscTryTypeMethod(*linesearch, destroy);
675648c30bcSBarry Smith   PetscCall(PetscViewerDestroy(&(*linesearch)->monitor));
676f4f49eeaSPierre Jolivet   PetscCall(SNESLineSearchMonitorCancel(*linesearch));
6779566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(linesearch));
6783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
679bf7f4e0aSPeter Brune }
680bf7f4e0aSPeter Brune 
681f40b411bSPeter Brune /*@
682dcf2fd19SBarry Smith   SNESLineSearchSetDefaultMonitor - Turns on/off printing useful information and debugging output about the line search.
683bf7f4e0aSPeter Brune 
684c3339decSBarry Smith   Logically Collective
685f6dfbefdSBarry Smith 
686bf7f4e0aSPeter Brune   Input Parameters:
687dcf2fd19SBarry Smith + linesearch - the linesearch object
68820f4b53cSBarry Smith - viewer     - an `PETSCVIEWERASCII` `PetscViewer` or `NULL` to turn off monitor
689bf7f4e0aSPeter Brune 
690f6dfbefdSBarry Smith   Options Database Key:
691dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - enables the monitor
692bf7f4e0aSPeter Brune 
693bf7f4e0aSPeter Brune   Level: intermediate
694bf7f4e0aSPeter Brune 
695e4094ef1SJacob Faibussowitsch   Developer Notes:
696f6dfbefdSBarry Smith   This monitor is implemented differently than the other line search monitors that are set with
697f6dfbefdSBarry Smith   `SNESLineSearchMonitorSet()` since it is called in many locations of the line search routines to display aspects of the
698d12e167eSBarry Smith   line search that are not visible to the other monitors.
699dcf2fd19SBarry Smith 
700420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `PETSCVIEWERASCII`, `SNESGetLineSearch()`, `SNESLineSearchGetDefaultMonitor()`, `PetscViewer`, `SNESLineSearchSetMonitor()`,
701f6dfbefdSBarry Smith           `SNESLineSearchMonitorSetFromOptions()`
702bf7f4e0aSPeter Brune @*/
703d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetDefaultMonitor(SNESLineSearch linesearch, PetscViewer viewer)
704d71ae5a4SJacob Faibussowitsch {
705bf7f4e0aSPeter Brune   PetscFunctionBegin;
706648c30bcSBarry Smith   PetscCall(PetscViewerDestroy(&linesearch->monitor));
707dcf2fd19SBarry Smith   linesearch->monitor = viewer;
7083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
709bf7f4e0aSPeter Brune }
710bf7f4e0aSPeter Brune 
711f40b411bSPeter Brune /*@
712420bcc1bSBarry Smith   SNESLineSearchGetDefaultMonitor - Gets the `PetscViewer` instance for the default line search monitor that is turned on with `SNESLineSearchSetDefaultMonitor()`
713f6dfbefdSBarry Smith 
714c3339decSBarry Smith   Logically Collective
7156a388c36SPeter Brune 
716f190f2fcSBarry Smith   Input Parameter:
717420bcc1bSBarry Smith . linesearch - the line search context
718f40b411bSPeter Brune 
719f190f2fcSBarry Smith   Output Parameter:
7208e557f58SPeter Brune . monitor - monitor context
721f40b411bSPeter Brune 
722f40b411bSPeter Brune   Level: intermediate
723f40b411bSPeter Brune 
724420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchSetDefaultMonitor()`, `PetscViewer`
725f40b411bSPeter Brune @*/
726d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetDefaultMonitor(SNESLineSearch linesearch, PetscViewer *monitor)
727d71ae5a4SJacob Faibussowitsch {
7286a388c36SPeter Brune   PetscFunctionBegin;
729f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
7306a388c36SPeter Brune   *monitor = linesearch->monitor;
7313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7326a388c36SPeter Brune }
7336a388c36SPeter Brune 
734dcf2fd19SBarry Smith /*@C
735420bcc1bSBarry Smith   SNESLineSearchMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated in the options database
736dcf2fd19SBarry Smith 
737c3339decSBarry Smith   Collective
738dcf2fd19SBarry Smith 
739dcf2fd19SBarry Smith   Input Parameters:
740420bcc1bSBarry Smith + ls           - `SNESLineSearch` object to monitor
741420bcc1bSBarry Smith . name         - the monitor type
742dcf2fd19SBarry Smith . help         - message indicating what monitoring is done
743dcf2fd19SBarry Smith . manual       - manual page for the monitor
744*49abdd8aSBarry Smith . monitor      - the monitor function, must use `PetscViewerAndFormat` as its context
745f6dfbefdSBarry Smith - monitorsetup - a function that is called once ONLY if the user selected this monitor that may set additional features of the `SNESLineSearch` or `PetscViewer`
746dcf2fd19SBarry Smith 
747420bcc1bSBarry Smith   Calling sequence of `monitor`:
748420bcc1bSBarry Smith + ls - `SNESLineSearch` object being monitored
749420bcc1bSBarry Smith - vf - a `PetscViewerAndFormat` struct that provides the `PetscViewer` and `PetscViewerFormat` being used
750dcf2fd19SBarry Smith 
751420bcc1bSBarry Smith   Calling sequence of `monitorsetup`:
752420bcc1bSBarry Smith + ls - `SNESLineSearch` object being monitored
753420bcc1bSBarry Smith - vf - a `PetscViewerAndFormat` struct that provides the `PetscViewer` and `PetscViewerFormat` being used
754420bcc1bSBarry Smith 
755420bcc1bSBarry Smith   Level: advanced
756420bcc1bSBarry Smith 
757648c30bcSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetMonitor()`, `PetscOptionsCreateViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
758db781477SPatrick Sanan           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
759e4094ef1SJacob Faibussowitsch           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
760db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
761c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
762db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
763db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
764dcf2fd19SBarry Smith @*/
765420bcc1bSBarry Smith 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))
766d71ae5a4SJacob Faibussowitsch {
767dcf2fd19SBarry Smith   PetscViewer       viewer;
768dcf2fd19SBarry Smith   PetscViewerFormat format;
769dcf2fd19SBarry Smith   PetscBool         flg;
770dcf2fd19SBarry Smith 
771dcf2fd19SBarry Smith   PetscFunctionBegin;
772648c30bcSBarry Smith   PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)ls), ((PetscObject)ls)->options, ((PetscObject)ls)->prefix, name, &viewer, &format, &flg));
773dcf2fd19SBarry Smith   if (flg) {
774d12e167eSBarry Smith     PetscViewerAndFormat *vf;
7759566063dSJacob Faibussowitsch     PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf));
776648c30bcSBarry Smith     PetscCall(PetscViewerDestroy(&viewer));
7771baa6e33SBarry Smith     if (monitorsetup) PetscCall((*monitorsetup)(ls, vf));
778*49abdd8aSBarry Smith     PetscCall(SNESLineSearchMonitorSet(ls, (PetscErrorCode (*)(SNESLineSearch, void *))monitor, vf, (PetscCtxDestroyFn *)PetscViewerAndFormatDestroy));
779dcf2fd19SBarry Smith   }
7803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
781dcf2fd19SBarry Smith }
782dcf2fd19SBarry Smith 
783f40b411bSPeter Brune /*@
784f1c6b773SPeter Brune   SNESLineSearchSetFromOptions - Sets options for the line search
785f40b411bSPeter Brune 
786c3339decSBarry Smith   Logically Collective
787f6dfbefdSBarry Smith 
788f899ff85SJose E. Roman   Input Parameter:
789420bcc1bSBarry Smith . linesearch - a `SNESLineSearch` line search context
790f40b411bSPeter Brune 
791f40b411bSPeter Brune   Options Database Keys:
7920b00b554SBarry Smith + -snes_linesearch_type <type>                                      - basic (or equivalently none), bt, l2, cp, nleqerr, shell
7933c7d6663SPeter Brune . -snes_linesearch_order <order>                                    - 1, 2, 3.  Most types only support certain orders (bt supports 2 or 3)
794f6dfbefdSBarry Smith . -snes_linesearch_norms                                            - Turn on/off the linesearch norms for the basic linesearch typem (`SNESLineSearchSetComputeNorms()`)
79571eef1aeSPeter Brune . -snes_linesearch_minlambda                                        - The minimum step length
7961a9b3a06SPeter Brune . -snes_linesearch_maxstep                                          - The maximum step size
7971a9b3a06SPeter Brune . -snes_linesearch_rtol                                             - Relative tolerance for iterative line searches
7981a9b3a06SPeter Brune . -snes_linesearch_atol                                             - Absolute tolerance for iterative line searches
7991a9b3a06SPeter Brune . -snes_linesearch_ltol                                             - Change in lambda tolerance for iterative line searches
8001a9b3a06SPeter Brune . -snes_linesearch_max_it                                           - The number of iterations for iterative line searches
801dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename]                              - Print progress of line searches
802dcf2fd19SBarry Smith . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine
8038e557f58SPeter Brune . -snes_linesearch_damping                                          - The linesearch damping parameter
804cd7522eaSPeter Brune . -snes_linesearch_keeplambda                                       - Keep the previous search length as the initial guess.
8051a9b3a06SPeter Brune . -snes_linesearch_precheck_picard                                  - Use precheck that speeds up convergence of picard method
806d8d34be6SBarry Smith - -snes_linesearch_precheck_picard_angle                            - Angle used in Picard precheck method
807f40b411bSPeter Brune 
808f40b411bSPeter Brune   Level: intermediate
809f40b411bSPeter Brune 
810420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchCreate()`, `SNESLineSearchSetOrder()`, `SNESLineSearchSetType()`, `SNESLineSearchSetTolerances()`, `SNESLineSearchSetDamping()`, `SNESLineSearchPreCheckPicard()`,
811db781477SPatrick Sanan           `SNESLineSearchType`, `SNESLineSearchSetComputeNorms()`
812f40b411bSPeter Brune @*/
813d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch)
814d71ae5a4SJacob Faibussowitsch {
8151a4f838cSPeter Brune   const char *deft = SNESLINESEARCHBASIC;
816bf7f4e0aSPeter Brune   char        type[256];
817bf7f4e0aSPeter Brune   PetscBool   flg, set;
818dcf2fd19SBarry Smith   PetscViewer viewer;
819bf388a1fSBarry Smith 
820bf7f4e0aSPeter Brune   PetscFunctionBegin;
8219566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchRegisterAll());
822bf7f4e0aSPeter Brune 
823d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)linesearch);
824f5af7f23SKarl Rupp   if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name;
8259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-snes_linesearch_type", "Linesearch type", "SNESLineSearchSetType", SNESLineSearchList, deft, type, 256, &flg));
826bf7f4e0aSPeter Brune   if (flg) {
8279566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetType(linesearch, type));
828bf7f4e0aSPeter Brune   } else if (!((PetscObject)linesearch)->type_name) {
8299566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetType(linesearch, deft));
830bf7f4e0aSPeter Brune   }
831bf7f4e0aSPeter Brune 
832648c30bcSBarry Smith   PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)linesearch), ((PetscObject)linesearch)->options, ((PetscObject)linesearch)->prefix, "-snes_linesearch_monitor", &viewer, NULL, &set));
833648c30bcSBarry Smith   if (set) PetscCall(SNESLineSearchSetDefaultMonitor(linesearch, viewer));
8349566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchMonitorSetFromOptions(linesearch, "-snes_linesearch_monitor_solution_update", "View correction at each iteration", "SNESLineSearchMonitorSolutionUpdate", SNESLineSearchMonitorSolutionUpdate, NULL));
835bf7f4e0aSPeter Brune 
8361a9b3a06SPeter Brune   /* tolerances */
8379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_minlambda", "Minimum step length", "SNESLineSearchSetTolerances", linesearch->steptol, &linesearch->steptol, NULL));
8389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_maxstep", "Maximum step size", "SNESLineSearchSetTolerances", linesearch->maxstep, &linesearch->maxstep, NULL));
8399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_rtol", "Relative tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->rtol, &linesearch->rtol, NULL));
8409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_atol", "Absolute tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->atol, &linesearch->atol, NULL));
8419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_ltol", "Change in lambda tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->ltol, &linesearch->ltol, NULL));
8429566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_linesearch_max_it", "Maximum iterations for iterative line searches", "SNESLineSearchSetTolerances", linesearch->max_its, &linesearch->max_its, NULL));
8437a35526eSPeter Brune 
8441a9b3a06SPeter Brune   /* damping parameters */
8459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_damping", "Line search damping and initial step guess", "SNESLineSearchSetDamping", linesearch->damping, &linesearch->damping, NULL));
8461a9b3a06SPeter Brune 
8479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_linesearch_keeplambda", "Use previous lambda as damping", "SNESLineSearchSetKeepLambda", linesearch->keeplambda, &linesearch->keeplambda, NULL));
8481a9b3a06SPeter Brune 
8491a9b3a06SPeter Brune   /* precheck */
8509566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_linesearch_precheck_picard", "Use a correction that sometimes improves convergence of Picard iteration", "SNESLineSearchPreCheckPicard", flg, &flg, &set));
8517a35526eSPeter Brune   if (set) {
8527a35526eSPeter Brune     if (flg) {
8537a35526eSPeter Brune       linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */
854f5af7f23SKarl Rupp 
855d0609cedSBarry Smith       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));
8569566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchSetPreCheck(linesearch, SNESLineSearchPreCheckPicard, &linesearch->precheck_picard_angle));
8577a35526eSPeter Brune     } else {
8589566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchSetPreCheck(linesearch, NULL, NULL));
8597a35526eSPeter Brune     }
8607a35526eSPeter Brune   }
8619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_linesearch_order", "Order of approximation used in the line search", "SNESLineSearchSetOrder", linesearch->order, &linesearch->order, NULL));
8629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_linesearch_norms", "Compute final norms in line search", "SNESLineSearchSetComputeNorms", linesearch->norms, &linesearch->norms, NULL));
8637a35526eSPeter Brune 
864dbbe0bcdSBarry Smith   PetscTryTypeMethod(linesearch, setfromoptions, PetscOptionsObject);
8655a9b6599SPeter Brune 
866dbbe0bcdSBarry Smith   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)linesearch, PetscOptionsObject));
867d0609cedSBarry Smith   PetscOptionsEnd();
8683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
869bf7f4e0aSPeter Brune }
870bf7f4e0aSPeter Brune 
871f40b411bSPeter Brune /*@
872f190f2fcSBarry Smith   SNESLineSearchView - Prints useful information about the line search
873f40b411bSPeter Brune 
87420f4b53cSBarry Smith   Logically Collective
87520f4b53cSBarry Smith 
876f40b411bSPeter Brune   Input Parameters:
8772fe279fdSBarry Smith + linesearch - line search context
878420bcc1bSBarry Smith - viewer     - the `PetscViewer` to display the line search information to
879f40b411bSPeter Brune 
880f40b411bSPeter Brune   Level: intermediate
881f40b411bSPeter Brune 
882420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `PetscViewer`, `SNESLineSearchCreate()`
883f40b411bSPeter Brune @*/
884d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer)
885d71ae5a4SJacob Faibussowitsch {
8867f1410a3SPeter Brune   PetscBool iascii;
887bf388a1fSBarry Smith 
888bf7f4e0aSPeter Brune   PetscFunctionBegin;
8897f1410a3SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
89048a46eb9SPierre Jolivet   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch), &viewer));
8917f1410a3SPeter Brune   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
8927f1410a3SPeter Brune   PetscCheckSameComm(linesearch, 1, viewer, 2);
893f40b411bSPeter Brune 
8949566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
8957f1410a3SPeter Brune   if (iascii) {
8969566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)linesearch, viewer));
8979566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
898dbbe0bcdSBarry Smith     PetscTryTypeMethod(linesearch, view, viewer);
8999566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
9009566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep, (double)linesearch->steptol));
9019566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol, (double)linesearch->atol, (double)linesearch->ltol));
90263a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  maximum iterations=%" PetscInt_FMT "\n", linesearch->max_its));
9036b2b7091SBarry Smith     if (linesearch->ops->precheck) {
9046b2b7091SBarry Smith       if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) {
90563a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  using precheck step to speed up Picard convergence\n"));
9067f1410a3SPeter Brune       } else {
90763a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  using user-defined precheck step\n"));
9087f1410a3SPeter Brune       }
9097f1410a3SPeter Brune     }
91048a46eb9SPierre Jolivet     if (linesearch->ops->postcheck) PetscCall(PetscViewerASCIIPrintf(viewer, "  using user-defined postcheck step\n"));
9117f1410a3SPeter Brune   }
9123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
913bf7f4e0aSPeter Brune }
914bf7f4e0aSPeter Brune 
915cc4c1da9SBarry Smith /*@
916420bcc1bSBarry Smith   SNESLineSearchGetType - Gets the `SNESLinesearchType` of a `SNESLineSearch`
917a80ff896SJed Brown 
918c3339decSBarry Smith   Logically Collective
919a80ff896SJed Brown 
9202fe279fdSBarry Smith   Input Parameter:
921420bcc1bSBarry Smith . linesearch - the line search context
922a80ff896SJed Brown 
9232fe279fdSBarry Smith   Output Parameter:
9242fe279fdSBarry Smith . type - The type of line search, or `NULL` if not set
925a80ff896SJed Brown 
926a80ff896SJed Brown   Level: intermediate
927a80ff896SJed Brown 
928420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchSetType()`
929a80ff896SJed Brown @*/
930d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetType(SNESLineSearch linesearch, SNESLineSearchType *type)
931d71ae5a4SJacob Faibussowitsch {
932a80ff896SJed Brown   PetscFunctionBegin;
933a80ff896SJed Brown   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
9344f572ea9SToby Isaac   PetscAssertPointer(type, 2);
935a80ff896SJed Brown   *type = ((PetscObject)linesearch)->type_name;
9363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
937a80ff896SJed Brown }
938a80ff896SJed Brown 
939cc4c1da9SBarry Smith /*@
940420bcc1bSBarry Smith   SNESLineSearchSetType - Sets the `SNESLinesearchType` of a `SNESLineSearch`
941f40b411bSPeter Brune 
942c3339decSBarry Smith   Logically Collective
943f190f2fcSBarry Smith 
944f40b411bSPeter Brune   Input Parameters:
945420bcc1bSBarry Smith + linesearch - the line search context
946ceaaa498SBarry Smith - type       - The type of line search to be used, see `SNESLineSearchType`
9471fe24845SBarry Smith 
9483c7db156SBarry Smith   Options Database Key:
9490b00b554SBarry Smith . -snes_linesearch_type <type> - basic (or equivalently none), bt, l2, cp, nleqerr, shell
950cd7522eaSPeter Brune 
951f40b411bSPeter Brune   Level: intermediate
952f40b411bSPeter Brune 
953420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchGetType()`
954f40b411bSPeter Brune @*/
955d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type)
956d71ae5a4SJacob Faibussowitsch {
957bf7f4e0aSPeter Brune   PetscBool match;
9585f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(SNESLineSearch);
959bf7f4e0aSPeter Brune 
960bf7f4e0aSPeter Brune   PetscFunctionBegin;
961f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
9624f572ea9SToby Isaac   PetscAssertPointer(type, 2);
963bf7f4e0aSPeter Brune 
9649566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)linesearch, type, &match));
9653ba16761SJacob Faibussowitsch   if (match) PetscFunctionReturn(PETSC_SUCCESS);
966bf7f4e0aSPeter Brune 
9679566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(SNESLineSearchList, type, &r));
9686adde796SStefano Zampini   PetscCheck(r, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested Line Search type %s", type);
969bf7f4e0aSPeter Brune   /* Destroy the previous private line search context */
970213acdd3SPierre Jolivet   PetscTryTypeMethod(linesearch, destroy);
9710298fd71SBarry Smith   linesearch->ops->destroy = NULL;
972f1c6b773SPeter Brune   /* Reinitialize function pointers in SNESLineSearchOps structure */
9739e5d0892SLisandro Dalcin   linesearch->ops->apply          = NULL;
9749e5d0892SLisandro Dalcin   linesearch->ops->view           = NULL;
9759e5d0892SLisandro Dalcin   linesearch->ops->setfromoptions = NULL;
9769e5d0892SLisandro Dalcin   linesearch->ops->destroy        = NULL;
977bf7f4e0aSPeter Brune 
9789566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)linesearch, type));
9799566063dSJacob Faibussowitsch   PetscCall((*r)(linesearch));
9803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
981bf7f4e0aSPeter Brune }
982bf7f4e0aSPeter Brune 
983f40b411bSPeter Brune /*@
984f6dfbefdSBarry Smith   SNESLineSearchSetSNES - Sets the `SNES` for the linesearch for function evaluation.
985f40b411bSPeter Brune 
986f40b411bSPeter Brune   Input Parameters:
987420bcc1bSBarry Smith + linesearch - the line search context
988420bcc1bSBarry Smith - snes       - The `SNES` instance
989f40b411bSPeter Brune 
99078bcb3b5SPeter Brune   Level: developer
99178bcb3b5SPeter Brune 
992f6dfbefdSBarry Smith   Note:
993f190f2fcSBarry Smith   This happens automatically when the line search is obtained/created with
994f6dfbefdSBarry Smith   `SNESGetLineSearch()`.  This routine is therefore mainly called within `SNES`
99578bcb3b5SPeter Brune   implementations.
996f40b411bSPeter Brune 
997420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetSNES()`, `SNESLineSearchSetVecs()`
998f40b411bSPeter Brune @*/
999d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes)
1000d71ae5a4SJacob Faibussowitsch {
1001bf7f4e0aSPeter Brune   PetscFunctionBegin;
1002f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1003bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 2);
1004bf7f4e0aSPeter Brune   linesearch->snes = snes;
10053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1006bf7f4e0aSPeter Brune }
1007bf7f4e0aSPeter Brune 
1008f40b411bSPeter Brune /*@
1009f6dfbefdSBarry Smith   SNESLineSearchGetSNES - Gets the `SNES` instance associated with the line search.
1010f6dfbefdSBarry Smith 
1011f6dfbefdSBarry Smith   Not Collective
1012f40b411bSPeter Brune 
10132fe279fdSBarry Smith   Input Parameter:
1014420bcc1bSBarry Smith . linesearch - the line search context
1015f40b411bSPeter Brune 
10162fe279fdSBarry Smith   Output Parameter:
10172fe279fdSBarry Smith . snes - The `SNES` instance
1018f40b411bSPeter Brune 
10198141a3b9SPeter Brune   Level: developer
1020f40b411bSPeter Brune 
1021420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESType`, `SNESLineSearchSetVecs()`
1022f40b411bSPeter Brune @*/
1023d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes)
1024d71ae5a4SJacob Faibussowitsch {
1025bf7f4e0aSPeter Brune   PetscFunctionBegin;
1026f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
10274f572ea9SToby Isaac   PetscAssertPointer(snes, 2);
1028bf7f4e0aSPeter Brune   *snes = linesearch->snes;
10293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1030bf7f4e0aSPeter Brune }
1031bf7f4e0aSPeter Brune 
1032f40b411bSPeter Brune /*@
1033420bcc1bSBarry Smith   SNESLineSearchGetLambda - Gets the last line search steplength used
1034f40b411bSPeter Brune 
1035f6dfbefdSBarry Smith   Not Collective
1036f6dfbefdSBarry Smith 
10372fe279fdSBarry Smith   Input Parameter:
1038420bcc1bSBarry Smith . linesearch - the line search context
1039f40b411bSPeter Brune 
10402fe279fdSBarry Smith   Output Parameter:
1041f6dfbefdSBarry Smith . lambda - The last steplength computed during `SNESLineSearchApply()`
1042f40b411bSPeter Brune 
104378bcb3b5SPeter Brune   Level: advanced
104478bcb3b5SPeter Brune 
1045f6dfbefdSBarry Smith   Note:
10468e557f58SPeter Brune   This is useful in methods where the solver is ill-scaled and
104778bcb3b5SPeter Brune   requires some adaptive notion of the difference in scale between the
1048f6dfbefdSBarry Smith   solution and the function.  For instance, `SNESQN` may be scaled by the
104978bcb3b5SPeter Brune   line search lambda using the argument -snes_qn_scaling ls.
105078bcb3b5SPeter Brune 
1051420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetLambda()`, `SNESLineSearchGetDamping()`, `SNESLineSearchApply()`
1052f40b411bSPeter Brune @*/
1053d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch linesearch, PetscReal *lambda)
1054d71ae5a4SJacob Faibussowitsch {
10556a388c36SPeter Brune   PetscFunctionBegin;
1056f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
10574f572ea9SToby Isaac   PetscAssertPointer(lambda, 2);
10586a388c36SPeter Brune   *lambda = linesearch->lambda;
10593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10606a388c36SPeter Brune }
10616a388c36SPeter Brune 
1062f40b411bSPeter Brune /*@
1063f6dfbefdSBarry Smith   SNESLineSearchSetLambda - Sets the line search steplength
1064f40b411bSPeter Brune 
1065f40b411bSPeter Brune   Input Parameters:
10668e557f58SPeter Brune + linesearch - line search context
1067420bcc1bSBarry Smith - lambda     - The steplength to use
1068f40b411bSPeter Brune 
106920f4b53cSBarry Smith   Level: advanced
107020f4b53cSBarry Smith 
1071f6dfbefdSBarry Smith   Note:
1072f6dfbefdSBarry Smith   This routine is typically used within implementations of `SNESLineSearchApply()`
1073f6dfbefdSBarry Smith   to set the final steplength.  This routine (and `SNESLineSearchGetLambda()`) were
1074cd7522eaSPeter Brune   added in order to facilitate Quasi-Newton methods that use the previous steplength
1075cd7522eaSPeter Brune   as an inner scaling parameter.
1076cd7522eaSPeter Brune 
1077420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetLambda()`
1078f40b411bSPeter Brune @*/
1079d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
1080d71ae5a4SJacob Faibussowitsch {
10816a388c36SPeter Brune   PetscFunctionBegin;
1082f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
10836a388c36SPeter Brune   linesearch->lambda = lambda;
10843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10856a388c36SPeter Brune }
10866a388c36SPeter Brune 
1087f40b411bSPeter Brune /*@
1088ceaaa498SBarry Smith   SNESLineSearchGetTolerances - Gets the tolerances for the line search.
1089f40b411bSPeter Brune 
1090f6dfbefdSBarry Smith   Not Collective
1091f6dfbefdSBarry Smith 
1092f899ff85SJose E. Roman   Input Parameter:
1093420bcc1bSBarry Smith . linesearch - the line search context
1094f40b411bSPeter Brune 
1095f40b411bSPeter Brune   Output Parameters:
1096b13b64c2SBarry Smith + steptol - The minimum steplength
10976cc8e53bSPeter Brune . maxstep - The maximum steplength
1098516fe3c3SPeter Brune . rtol    - The relative tolerance for iterative line searches
1099516fe3c3SPeter Brune . atol    - The absolute tolerance for iterative line searches
1100516fe3c3SPeter Brune . ltol    - The change in lambda tolerance for iterative line searches
1101e4094ef1SJacob Faibussowitsch - max_its - The maximum number of iterations of the line search
1102f40b411bSPeter Brune 
110378bcb3b5SPeter Brune   Level: intermediate
110478bcb3b5SPeter Brune 
1105f6dfbefdSBarry Smith   Note:
110678bcb3b5SPeter Brune   Different line searches may implement these parameters slightly differently as
11073c7d6663SPeter Brune   the type requires.
1108516fe3c3SPeter Brune 
1109420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetTolerances()`
1110f40b411bSPeter Brune @*/
1111d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch linesearch, PetscReal *steptol, PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
1112d71ae5a4SJacob Faibussowitsch {
11136a388c36SPeter Brune   PetscFunctionBegin;
1114f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1115516fe3c3SPeter Brune   if (steptol) {
11164f572ea9SToby Isaac     PetscAssertPointer(steptol, 2);
11176a388c36SPeter Brune     *steptol = linesearch->steptol;
1118516fe3c3SPeter Brune   }
1119516fe3c3SPeter Brune   if (maxstep) {
11204f572ea9SToby Isaac     PetscAssertPointer(maxstep, 3);
1121b13b64c2SBarry Smith     *maxstep = linesearch->maxstep;
1122516fe3c3SPeter Brune   }
1123516fe3c3SPeter Brune   if (rtol) {
11244f572ea9SToby Isaac     PetscAssertPointer(rtol, 4);
1125516fe3c3SPeter Brune     *rtol = linesearch->rtol;
1126516fe3c3SPeter Brune   }
1127516fe3c3SPeter Brune   if (atol) {
11284f572ea9SToby Isaac     PetscAssertPointer(atol, 5);
1129516fe3c3SPeter Brune     *atol = linesearch->atol;
1130516fe3c3SPeter Brune   }
1131516fe3c3SPeter Brune   if (ltol) {
11324f572ea9SToby Isaac     PetscAssertPointer(ltol, 6);
1133516fe3c3SPeter Brune     *ltol = linesearch->ltol;
1134516fe3c3SPeter Brune   }
1135516fe3c3SPeter Brune   if (max_its) {
11364f572ea9SToby Isaac     PetscAssertPointer(max_its, 7);
1137516fe3c3SPeter Brune     *max_its = linesearch->max_its;
1138516fe3c3SPeter Brune   }
11393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11406a388c36SPeter Brune }
11416a388c36SPeter Brune 
1142f40b411bSPeter Brune /*@
1143ceaaa498SBarry Smith   SNESLineSearchSetTolerances -  Sets the tolerances for the linesearch.
1144f40b411bSPeter Brune 
1145f6dfbefdSBarry Smith   Collective
1146f6dfbefdSBarry Smith 
1147f40b411bSPeter Brune   Input Parameters:
1148420bcc1bSBarry Smith + linesearch - the line search context
1149516fe3c3SPeter Brune . steptol    - The minimum steplength
11506cc8e53bSPeter Brune . maxstep    - The maximum steplength
1151516fe3c3SPeter Brune . rtol       - The relative tolerance for iterative line searches
1152516fe3c3SPeter Brune . atol       - The absolute tolerance for iterative line searches
1153516fe3c3SPeter Brune . ltol       - The change in lambda tolerance for iterative line searches
1154420bcc1bSBarry Smith - max_it     - The maximum number of iterations of the line search
1155420bcc1bSBarry Smith 
1156420bcc1bSBarry Smith   Options Database Keys:
1157420bcc1bSBarry Smith + -snes_linesearch_minlambda - The minimum step length
1158420bcc1bSBarry Smith . -snes_linesearch_maxstep   - The maximum step size
1159420bcc1bSBarry Smith . -snes_linesearch_rtol      - Relative tolerance for iterative line searches
1160420bcc1bSBarry Smith . -snes_linesearch_atol      - Absolute tolerance for iterative line searches
1161420bcc1bSBarry Smith . -snes_linesearch_ltol      - Change in lambda tolerance for iterative line searches
1162420bcc1bSBarry Smith - -snes_linesearch_max_it    - The number of iterations for iterative line searches
1163f40b411bSPeter Brune 
116420f4b53cSBarry Smith   Level: intermediate
116520f4b53cSBarry Smith 
1166f6dfbefdSBarry Smith   Note:
1167420bcc1bSBarry Smith   The user may choose to not set any of the tolerances using `PETSC_DEFAULT` in place of an argument.
1168f40b411bSPeter Brune 
1169420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetTolerances()`
1170f40b411bSPeter Brune @*/
1171420bcc1bSBarry Smith PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch linesearch, PetscReal steptol, PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_it)
1172d71ae5a4SJacob Faibussowitsch {
11736a388c36SPeter Brune   PetscFunctionBegin;
1174f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1175d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, steptol, 2);
1176d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, maxstep, 3);
1177d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, rtol, 4);
1178d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, atol, 5);
1179d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, ltol, 6);
1180420bcc1bSBarry Smith   PetscValidLogicalCollectiveInt(linesearch, max_it, 7);
1181d3952184SSatish Balay 
118213bcc0bdSJacob Faibussowitsch   if (steptol != (PetscReal)PETSC_DEFAULT) {
11835f80ce2aSJacob Faibussowitsch     PetscCheck(steptol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Minimum step length %14.12e must be non-negative", (double)steptol);
11846a388c36SPeter Brune     linesearch->steptol = steptol;
1185d3952184SSatish Balay   }
1186d3952184SSatish Balay 
118713bcc0bdSJacob Faibussowitsch   if (maxstep != (PetscReal)PETSC_DEFAULT) {
11885f80ce2aSJacob Faibussowitsch     PetscCheck(maxstep >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum step length %14.12e must be non-negative", (double)maxstep);
1189516fe3c3SPeter Brune     linesearch->maxstep = maxstep;
1190d3952184SSatish Balay   }
1191d3952184SSatish Balay 
119213bcc0bdSJacob Faibussowitsch   if (rtol != (PetscReal)PETSC_DEFAULT) {
11932061ca32SJunchao Zhang     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);
1194516fe3c3SPeter Brune     linesearch->rtol = rtol;
1195d3952184SSatish Balay   }
1196d3952184SSatish Balay 
119713bcc0bdSJacob Faibussowitsch   if (atol != (PetscReal)PETSC_DEFAULT) {
11985f80ce2aSJacob Faibussowitsch     PetscCheck(atol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Absolute tolerance %14.12e must be non-negative", (double)atol);
1199516fe3c3SPeter Brune     linesearch->atol = atol;
1200d3952184SSatish Balay   }
1201d3952184SSatish Balay 
120213bcc0bdSJacob Faibussowitsch   if (ltol != (PetscReal)PETSC_DEFAULT) {
12035f80ce2aSJacob Faibussowitsch     PetscCheck(ltol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Lambda tolerance %14.12e must be non-negative", (double)ltol);
1204516fe3c3SPeter Brune     linesearch->ltol = ltol;
1205d3952184SSatish Balay   }
1206d3952184SSatish Balay 
1207420bcc1bSBarry Smith   if (max_it != PETSC_DEFAULT) {
1208420bcc1bSBarry Smith     PetscCheck(max_it >= 0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of iterations %" PetscInt_FMT " must be non-negative", max_it);
1209420bcc1bSBarry Smith     linesearch->max_its = max_it;
1210d3952184SSatish Balay   }
12113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12126a388c36SPeter Brune }
12136a388c36SPeter Brune 
1214f40b411bSPeter Brune /*@
1215f1c6b773SPeter Brune   SNESLineSearchGetDamping - Gets the line search damping parameter.
1216f40b411bSPeter Brune 
12172fe279fdSBarry Smith   Input Parameter:
1218420bcc1bSBarry Smith . linesearch - the line search context
1219f40b411bSPeter Brune 
12202fe279fdSBarry Smith   Output Parameter:
12218e557f58SPeter Brune . damping - The damping parameter
1222f40b411bSPeter Brune 
122378bcb3b5SPeter Brune   Level: advanced
1224f40b411bSPeter Brune 
1225420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearchGetStepTolerance()`, `SNESQN`
1226f40b411bSPeter Brune @*/
1227d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch linesearch, PetscReal *damping)
1228d71ae5a4SJacob Faibussowitsch {
12296a388c36SPeter Brune   PetscFunctionBegin;
1230f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
12314f572ea9SToby Isaac   PetscAssertPointer(damping, 2);
12326a388c36SPeter Brune   *damping = linesearch->damping;
12333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12346a388c36SPeter Brune }
12356a388c36SPeter Brune 
1236f40b411bSPeter Brune /*@
1237fd292e60Sprj-   SNESLineSearchSetDamping - Sets the line search damping parameter.
1238f40b411bSPeter Brune 
1239f40b411bSPeter Brune   Input Parameters:
1240420bcc1bSBarry Smith + linesearch - the line search context
124103fd4120SBarry Smith - damping    - The damping parameter
1242f40b411bSPeter Brune 
1243f6dfbefdSBarry Smith   Options Database Key:
1244f6dfbefdSBarry Smith . -snes_linesearch_damping <damping> - the damping value
1245f6dfbefdSBarry Smith 
1246f40b411bSPeter Brune   Level: intermediate
1247f40b411bSPeter Brune 
1248f6dfbefdSBarry Smith   Note:
1249f6dfbefdSBarry Smith   The `SNESLINESEARCHNONE` line search merely takes the update step scaled by the damping parameter.
1250420bcc1bSBarry Smith   The use of the damping parameter in the `SNESLINESEARCHL2` and `SNESLINESEARCHCP` line searches is much more subtle;
125178bcb3b5SPeter Brune   it is used as a starting point in calculating the secant step. However, the eventual
1252420bcc1bSBarry Smith   step may be of greater length than the damping parameter.  In the `SNESLINESEARCHBT` line search it is
1253420bcc1bSBarry Smith   used as the maximum possible step length, as the `SNESLINESEARCHBT` line search only backtracks.
1254cd7522eaSPeter Brune 
1255420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetDamping()`
1256f40b411bSPeter Brune @*/
1257d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch linesearch, PetscReal damping)
1258d71ae5a4SJacob Faibussowitsch {
12596a388c36SPeter Brune   PetscFunctionBegin;
1260f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
12616a388c36SPeter Brune   linesearch->damping = damping;
12623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12636a388c36SPeter Brune }
12646a388c36SPeter Brune 
126559405d5eSPeter Brune /*@
126659405d5eSPeter Brune   SNESLineSearchGetOrder - Gets the line search approximation order.
126759405d5eSPeter Brune 
1268f6dfbefdSBarry Smith   Input Parameter:
1269420bcc1bSBarry Smith . linesearch - the line search context
127059405d5eSPeter Brune 
1271f6dfbefdSBarry Smith   Output Parameter:
12728e557f58SPeter Brune . order - The order
127359405d5eSPeter Brune 
127459405d5eSPeter Brune   Level: intermediate
127559405d5eSPeter Brune 
1276420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetOrder()`
127759405d5eSPeter Brune @*/
1278d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch linesearch, PetscInt *order)
1279d71ae5a4SJacob Faibussowitsch {
128059405d5eSPeter Brune   PetscFunctionBegin;
128159405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
12824f572ea9SToby Isaac   PetscAssertPointer(order, 2);
128359405d5eSPeter Brune   *order = linesearch->order;
12843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
128559405d5eSPeter Brune }
128659405d5eSPeter Brune 
128759405d5eSPeter Brune /*@
12881f8196a2SJed Brown   SNESLineSearchSetOrder - Sets the maximum order of the polynomial fit used in the line search
128959405d5eSPeter Brune 
129059405d5eSPeter Brune   Input Parameters:
1291420bcc1bSBarry Smith + linesearch - the line search context
1292ceaaa498SBarry Smith - order      - The order
129359405d5eSPeter Brune 
129459405d5eSPeter Brune   Level: intermediate
129559405d5eSPeter Brune 
1296ceaaa498SBarry Smith   Values for `order`\:
1297f6dfbefdSBarry Smith +  1 or `SNES_LINESEARCH_ORDER_LINEAR` - linear order
1298f6dfbefdSBarry Smith .  2 or `SNES_LINESEARCH_ORDER_QUADRATIC` - quadratic order
1299f6dfbefdSBarry Smith -  3 or `SNES_LINESEARCH_ORDER_CUBIC` - cubic order
130078bcb3b5SPeter Brune 
1301420bcc1bSBarry Smith   Options Database Key:
1302420bcc1bSBarry Smith . -snes_linesearch_order <order> - 1, 2, 3.  Most types only support certain orders (`SNESLINESEARCHBT` supports 2 or 3)
1303420bcc1bSBarry Smith 
1304ceaaa498SBarry Smith   Note:
1305ceaaa498SBarry Smith   These orders are supported by `SNESLINESEARCHBT` and `SNESLINESEARCHCP`
130659405d5eSPeter Brune 
1307420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetOrder()`, `SNESLineSearchSetDamping()`
130859405d5eSPeter Brune @*/
1309d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch linesearch, PetscInt order)
1310d71ae5a4SJacob Faibussowitsch {
131159405d5eSPeter Brune   PetscFunctionBegin;
131259405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
131359405d5eSPeter Brune   linesearch->order = order;
13143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
131559405d5eSPeter Brune }
131659405d5eSPeter Brune 
1317f40b411bSPeter Brune /*@
1318420bcc1bSBarry Smith   SNESLineSearchGetNorms - Gets the norms for the current solution `X`, the current update `Y`, and the current function value `F`.
1319f40b411bSPeter Brune 
1320f6dfbefdSBarry Smith   Not Collective
1321f6dfbefdSBarry Smith 
1322f899ff85SJose E. Roman   Input Parameter:
1323420bcc1bSBarry Smith . linesearch - the line search context
1324f40b411bSPeter Brune 
1325f40b411bSPeter Brune   Output Parameters:
1326f40b411bSPeter Brune + xnorm - The norm of the current solution
1327a68bbae5SBarry Smith . fnorm - The norm of the current function, this is the `norm(function(X))` where `X` is the current solution.
13286b095a85SStefano Zampini - ynorm - The norm of the current update (after scaling by the linesearch computed lambda)
1329f40b411bSPeter Brune 
133078bcb3b5SPeter Brune   Level: developer
1331f40b411bSPeter Brune 
1332420bcc1bSBarry Smith   Notes:
1333420bcc1bSBarry Smith   Some values may not be up-to-date at particular points in the code.
1334a68bbae5SBarry Smith 
1335a68bbae5SBarry Smith   This, in combination with `SNESLineSearchSetNorms()`, allow the line search and the `SNESSolve_XXX()` to share
1336a68bbae5SBarry Smith   computed values.
1337a68bbae5SBarry Smith 
13380241b274SPierre Jolivet .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetNorms()`, `SNESLineSearchGetVecs()`
1339f40b411bSPeter Brune @*/
1340d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal *xnorm, PetscReal *fnorm, PetscReal *ynorm)
1341d71ae5a4SJacob Faibussowitsch {
1342bf7f4e0aSPeter Brune   PetscFunctionBegin;
1343f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1344f5af7f23SKarl Rupp   if (xnorm) *xnorm = linesearch->xnorm;
1345f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
1346f5af7f23SKarl Rupp   if (ynorm) *ynorm = linesearch->ynorm;
13473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1348bf7f4e0aSPeter Brune }
1349bf7f4e0aSPeter Brune 
1350f40b411bSPeter Brune /*@
13511bd63e3eSSatish Balay   SNESLineSearchSetNorms - Sets the computed norms for the current solution `X`, the current update `Y`, and the current function value `F`.
1352f40b411bSPeter Brune 
1353c3339decSBarry Smith   Collective
1354f6dfbefdSBarry Smith 
1355f40b411bSPeter Brune   Input Parameters:
1356420bcc1bSBarry Smith + linesearch - the line search context
1357f40b411bSPeter Brune . xnorm      - The norm of the current solution
1358a68bbae5SBarry Smith . fnorm      - The norm of the current function, this is the `norm(function(X))` where `X` is the current solution
13596b095a85SStefano Zampini - ynorm      - The norm of the current update (after scaling by the linesearch computed lambda)
1360f40b411bSPeter Brune 
1361f6dfbefdSBarry Smith   Level: developer
1362f40b411bSPeter Brune 
1363420bcc1bSBarry Smith   Note:
1364420bcc1bSBarry Smith   This is called by the line search routines to store the values they have just computed
1365420bcc1bSBarry Smith 
1366420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()`
1367f40b411bSPeter Brune @*/
1368d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
1369d71ae5a4SJacob Faibussowitsch {
13706a388c36SPeter Brune   PetscFunctionBegin;
1371f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
13726a388c36SPeter Brune   linesearch->xnorm = xnorm;
13736a388c36SPeter Brune   linesearch->fnorm = fnorm;
13746a388c36SPeter Brune   linesearch->ynorm = ynorm;
13753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13766a388c36SPeter Brune }
13776a388c36SPeter Brune 
1378f40b411bSPeter Brune /*@
1379420bcc1bSBarry Smith   SNESLineSearchComputeNorms - Explicitly computes the norms of the current solution `X`, the current update `Y`, and the current function value `F`.
1380f40b411bSPeter Brune 
13812fe279fdSBarry Smith   Input Parameter:
1382420bcc1bSBarry Smith . linesearch - the line search context
1383f40b411bSPeter Brune 
138420f4b53cSBarry Smith   Options Database Key:
13858e557f58SPeter Brune . -snes_linesearch_norms - turn norm computation on or off
1386f40b411bSPeter Brune 
1387f40b411bSPeter Brune   Level: intermediate
1388f40b411bSPeter Brune 
1389420bcc1bSBarry Smith   Developer Note:
1390420bcc1bSBarry Smith   The options database key is misnamed. It should be -snes_linesearch_compute_norms
1391420bcc1bSBarry Smith 
1392420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms`, `SNESLineSearchSetNorms()`, `SNESLineSearchSetComputeNorms()`
1393f40b411bSPeter Brune @*/
1394d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
1395d71ae5a4SJacob Faibussowitsch {
13969bd66eb0SPeter Brune   SNES snes;
1397bf388a1fSBarry Smith 
13986a388c36SPeter Brune   PetscFunctionBegin;
13996a388c36SPeter Brune   if (linesearch->norms) {
14009bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
14019566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
14029566063dSJacob Faibussowitsch       PetscCall(VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
14039566063dSJacob Faibussowitsch       PetscCall(VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm));
14049566063dSJacob Faibussowitsch       PetscCall((*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm));
14059bd66eb0SPeter Brune     } else {
14069566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(linesearch->vec_func, NORM_2, &linesearch->fnorm));
14079566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
14089566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm));
14099566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(linesearch->vec_func, NORM_2, &linesearch->fnorm));
14109566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
14119566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(linesearch->vec_update, NORM_2, &linesearch->ynorm));
14126a388c36SPeter Brune     }
14139bd66eb0SPeter Brune   }
14143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14156a388c36SPeter Brune }
14166a388c36SPeter Brune 
14176f263ca3SPeter Brune /*@
14186f263ca3SPeter Brune   SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search.
14196f263ca3SPeter Brune 
14206f263ca3SPeter Brune   Input Parameters:
1421420bcc1bSBarry Smith + linesearch - the line search context
142278bcb3b5SPeter Brune - flg        - indicates whether or not to compute norms
14236f263ca3SPeter Brune 
142420f4b53cSBarry Smith   Options Database Key:
1425420bcc1bSBarry Smith . -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic (none) `SNESLINESEARCHBASIC` line search
14266f263ca3SPeter Brune 
142720f4b53cSBarry Smith   Level: intermediate
142820f4b53cSBarry Smith 
1429f6dfbefdSBarry Smith   Note:
1430f6dfbefdSBarry Smith   This is most relevant to the `SNESLINESEARCHBASIC` (or equivalently `SNESLINESEARCHNONE`) line search type since most line searches have a stopping criteria involving the norm.
14316f263ca3SPeter Brune 
1432420bcc1bSBarry Smith   Developer Note:
1433420bcc1bSBarry Smith   The options database key is misnamed. It should be -snes_linesearch_compute_norms
1434420bcc1bSBarry Smith 
1435420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetNorms()`, `SNESLineSearchComputeNorms()`, `SNESLINESEARCHBASIC`
14366f263ca3SPeter Brune @*/
1437d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg)
1438d71ae5a4SJacob Faibussowitsch {
14396f263ca3SPeter Brune   PetscFunctionBegin;
14406f263ca3SPeter Brune   linesearch->norms = flg;
14413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14426f263ca3SPeter Brune }
14436f263ca3SPeter Brune 
1444f40b411bSPeter Brune /*@
1445f6dfbefdSBarry Smith   SNESLineSearchGetVecs - Gets the vectors from the `SNESLineSearch` context
1446f6dfbefdSBarry Smith 
14478f14a041SBarry Smith   Not Collective but the vectors are parallel
1448f40b411bSPeter Brune 
1449f899ff85SJose E. Roman   Input Parameter:
1450420bcc1bSBarry Smith . linesearch - the line search context
1451f40b411bSPeter Brune 
1452f40b411bSPeter Brune   Output Parameters:
14536232e825SPeter Brune + X - Solution vector
14546232e825SPeter Brune . F - Function vector
14556232e825SPeter Brune . Y - Search direction vector
14566232e825SPeter Brune . W - Solution work vector
14576232e825SPeter Brune - G - Function work vector
14586232e825SPeter Brune 
145920f4b53cSBarry Smith   Level: advanced
146020f4b53cSBarry Smith 
14617bba9028SPeter Brune   Notes:
146220f4b53cSBarry Smith   At the beginning of a line search application, `X` should contain a
146320f4b53cSBarry Smith   solution and the vector `F` the function computed at `X`.  At the end of the
146420f4b53cSBarry Smith   line search application, `X` should contain the new solution, and `F` the
14656232e825SPeter Brune   function evaluated at the new solution.
1466f40b411bSPeter Brune 
1467f6dfbefdSBarry Smith   These vectors are owned by the `SNESLineSearch` and should not be destroyed by the caller
14682a7a6963SBarry Smith 
1469420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()`
1470f40b411bSPeter Brune @*/
1471d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch, Vec *X, Vec *F, Vec *Y, Vec *W, Vec *G)
1472d71ae5a4SJacob Faibussowitsch {
14736a388c36SPeter Brune   PetscFunctionBegin;
1474f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
14756a388c36SPeter Brune   if (X) {
14764f572ea9SToby Isaac     PetscAssertPointer(X, 2);
14776a388c36SPeter Brune     *X = linesearch->vec_sol;
14786a388c36SPeter Brune   }
14796a388c36SPeter Brune   if (F) {
14804f572ea9SToby Isaac     PetscAssertPointer(F, 3);
14816a388c36SPeter Brune     *F = linesearch->vec_func;
14826a388c36SPeter Brune   }
14836a388c36SPeter Brune   if (Y) {
14844f572ea9SToby Isaac     PetscAssertPointer(Y, 4);
14856a388c36SPeter Brune     *Y = linesearch->vec_update;
14866a388c36SPeter Brune   }
14876a388c36SPeter Brune   if (W) {
14884f572ea9SToby Isaac     PetscAssertPointer(W, 5);
14896a388c36SPeter Brune     *W = linesearch->vec_sol_new;
14906a388c36SPeter Brune   }
14916a388c36SPeter Brune   if (G) {
14924f572ea9SToby Isaac     PetscAssertPointer(G, 6);
14936a388c36SPeter Brune     *G = linesearch->vec_func_new;
14946a388c36SPeter Brune   }
14953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14966a388c36SPeter Brune }
14976a388c36SPeter Brune 
1498f40b411bSPeter Brune /*@
1499f6dfbefdSBarry Smith   SNESLineSearchSetVecs - Sets the vectors on the `SNESLineSearch` context
1500f6dfbefdSBarry Smith 
1501c3339decSBarry Smith   Logically Collective
1502f40b411bSPeter Brune 
1503f40b411bSPeter Brune   Input Parameters:
1504420bcc1bSBarry Smith + linesearch - the line search context
15056232e825SPeter Brune . X          - Solution vector
15066232e825SPeter Brune . F          - Function vector
15076232e825SPeter Brune . Y          - Search direction vector
15086232e825SPeter Brune . W          - Solution work vector
15096232e825SPeter Brune - G          - Function work vector
1510f40b411bSPeter Brune 
1511420bcc1bSBarry Smith   Level: developer
1512f40b411bSPeter Brune 
1513420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetNorms()`, `SNESLineSearchGetVecs()`
1514f40b411bSPeter Brune @*/
1515d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch, Vec X, Vec F, Vec Y, Vec W, Vec G)
1516d71ae5a4SJacob Faibussowitsch {
15176a388c36SPeter Brune   PetscFunctionBegin;
1518f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
15196a388c36SPeter Brune   if (X) {
15206a388c36SPeter Brune     PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
15216a388c36SPeter Brune     linesearch->vec_sol = X;
15226a388c36SPeter Brune   }
15236a388c36SPeter Brune   if (F) {
15246a388c36SPeter Brune     PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
15256a388c36SPeter Brune     linesearch->vec_func = F;
15266a388c36SPeter Brune   }
15276a388c36SPeter Brune   if (Y) {
15286a388c36SPeter Brune     PetscValidHeaderSpecific(Y, VEC_CLASSID, 4);
15296a388c36SPeter Brune     linesearch->vec_update = Y;
15306a388c36SPeter Brune   }
15316a388c36SPeter Brune   if (W) {
15326a388c36SPeter Brune     PetscValidHeaderSpecific(W, VEC_CLASSID, 5);
15336a388c36SPeter Brune     linesearch->vec_sol_new = W;
15346a388c36SPeter Brune   }
15356a388c36SPeter Brune   if (G) {
15366a388c36SPeter Brune     PetscValidHeaderSpecific(G, VEC_CLASSID, 6);
15376a388c36SPeter Brune     linesearch->vec_func_new = G;
15386a388c36SPeter Brune   }
15393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15406a388c36SPeter Brune }
15416a388c36SPeter Brune 
1542cc4c1da9SBarry Smith /*@
1543f1c6b773SPeter Brune   SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1544f6dfbefdSBarry Smith   `SNESLineSearch` options in the database.
1545e7058c64SPeter Brune 
1546c3339decSBarry Smith   Logically Collective
1547e7058c64SPeter Brune 
1548e7058c64SPeter Brune   Input Parameters:
1549f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context
1550e7058c64SPeter Brune - prefix     - the prefix to prepend to all option names
1551e7058c64SPeter Brune 
155220f4b53cSBarry Smith   Level: advanced
155320f4b53cSBarry Smith 
1554f6dfbefdSBarry Smith   Note:
1555e7058c64SPeter Brune   A hyphen (-) must NOT be given at the beginning of the prefix name.
1556e7058c64SPeter Brune   The first character of all runtime options is AUTOMATICALLY the hyphen.
1557e7058c64SPeter Brune 
1558420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch()`, `SNESLineSearchSetFromOptions()`, `SNESGetOptionsPrefix()`
1559e7058c64SPeter Brune @*/
1560d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch, const char prefix[])
1561d71ae5a4SJacob Faibussowitsch {
1562e7058c64SPeter Brune   PetscFunctionBegin;
1563f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
15649566063dSJacob Faibussowitsch   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)linesearch, prefix));
15653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1566e7058c64SPeter Brune }
1567e7058c64SPeter Brune 
1568cc4c1da9SBarry Smith /*@
1569f6dfbefdSBarry Smith   SNESLineSearchGetOptionsPrefix - Gets the prefix used for searching for all
1570f1c6b773SPeter Brune   SNESLineSearch options in the database.
1571e7058c64SPeter Brune 
1572e7058c64SPeter Brune   Not Collective
1573e7058c64SPeter Brune 
1574e7058c64SPeter Brune   Input Parameter:
1575f6dfbefdSBarry Smith . linesearch - the `SNESLineSearch` context
1576e7058c64SPeter Brune 
1577e7058c64SPeter Brune   Output Parameter:
1578e7058c64SPeter Brune . prefix - pointer to the prefix string used
1579e7058c64SPeter Brune 
1580e7058c64SPeter Brune   Level: advanced
1581e7058c64SPeter Brune 
1582e4094ef1SJacob Faibussowitsch   Fortran Notes:
158320f4b53cSBarry Smith   The user should pass in a string 'prefix' of
158420f4b53cSBarry Smith   sufficient length to hold the prefix.
158520f4b53cSBarry Smith 
1586420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESAppendOptionsPrefix()`
1587e7058c64SPeter Brune @*/
1588d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch, const char *prefix[])
1589d71ae5a4SJacob Faibussowitsch {
1590e7058c64SPeter Brune   PetscFunctionBegin;
1591f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
15929566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)linesearch, prefix));
15933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1594e7058c64SPeter Brune }
1595bf7f4e0aSPeter Brune 
15968d359177SBarry Smith /*@C
1597f6dfbefdSBarry Smith   SNESLineSearchSetWorkVecs - Sets work vectors for the line search.
1598f40b411bSPeter Brune 
1599d8d19677SJose E. Roman   Input Parameters:
1600f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context
1601f40b411bSPeter Brune - nwork      - the number of work vectors
1602f40b411bSPeter Brune 
1603f40b411bSPeter Brune   Level: developer
1604f40b411bSPeter Brune 
1605420bcc1bSBarry Smith   Developer Note:
1606420bcc1bSBarry Smith   This is called from within the set up routines for each of the line search types `SNESLineSearchType`
1607420bcc1bSBarry Smith 
1608420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESSetWorkVecs()`
1609f40b411bSPeter Brune @*/
1610d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork)
1611d71ae5a4SJacob Faibussowitsch {
1612bf7f4e0aSPeter Brune   PetscFunctionBegin;
16130fdf79fbSJacob Faibussowitsch   PetscCheck(linesearch->vec_sol, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
16149566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work));
16153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1616bf7f4e0aSPeter Brune }
1617bf7f4e0aSPeter Brune 
1618f40b411bSPeter Brune /*@
1619422a814eSBarry Smith   SNESLineSearchGetReason - Gets the success/failure status of the last line search application
1620f40b411bSPeter Brune 
16212fe279fdSBarry Smith   Input Parameter:
1622420bcc1bSBarry Smith . linesearch - the line search context
1623f40b411bSPeter Brune 
16242fe279fdSBarry Smith   Output Parameter:
1625422a814eSBarry Smith . result - The success or failure status
1626f40b411bSPeter Brune 
162720f4b53cSBarry Smith   Level: developer
162820f4b53cSBarry Smith 
1629f6dfbefdSBarry Smith   Note:
1630420bcc1bSBarry Smith   This is typically called after `SNESLineSearchApply()` in order to determine if the line search failed
1631420bcc1bSBarry Smith   (and set into the `SNES` convergence accordingly).
1632cd7522eaSPeter Brune 
1633420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetReason()`, `SNESLineSearchReason`
1634f40b411bSPeter Brune @*/
1635d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result)
1636d71ae5a4SJacob Faibussowitsch {
1637bf7f4e0aSPeter Brune   PetscFunctionBegin;
1638f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
16394f572ea9SToby Isaac   PetscAssertPointer(result, 2);
1640422a814eSBarry Smith   *result = linesearch->result;
16413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1642bf7f4e0aSPeter Brune }
1643bf7f4e0aSPeter Brune 
1644f40b411bSPeter Brune /*@
1645420bcc1bSBarry Smith   SNESLineSearchSetReason - Sets the success/failure status of the line search application
1646f40b411bSPeter Brune 
16475d83a8b1SBarry Smith   Logically Collective; No Fortran Support
16485d83a8b1SBarry Smith 
1649f40b411bSPeter Brune   Input Parameters:
1650420bcc1bSBarry Smith + linesearch - the line search context
1651422a814eSBarry Smith - result     - The success or failure status
1652f40b411bSPeter Brune 
165320f4b53cSBarry Smith   Level: developer
165420f4b53cSBarry Smith 
1655f6dfbefdSBarry Smith   Note:
1656420bcc1bSBarry Smith   This is typically called in a `SNESLineSearchType` implementation of `SNESLineSearchApply()` or a `SNESLINESEARCHSHELL` implementation to set
1657cd7522eaSPeter Brune   the success or failure of the line search method.
1658cd7522eaSPeter Brune 
1659420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchReason`, `SNESLineSearchGetSResult()`
1660f40b411bSPeter Brune @*/
1661d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result)
1662d71ae5a4SJacob Faibussowitsch {
16636a388c36SPeter Brune   PetscFunctionBegin;
16645fd66863SKarl Rupp   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1665422a814eSBarry Smith   linesearch->result = result;
16663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16676a388c36SPeter Brune }
16686a388c36SPeter Brune 
1669ceaaa498SBarry Smith // PetscClangLinter pragma disable: -fdoc-param-list-func-parameter-documentation
16709bd66eb0SPeter Brune /*@C
1671f1c6b773SPeter Brune   SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.
16729bd66eb0SPeter Brune 
1673c3339decSBarry Smith   Logically Collective
1674f6dfbefdSBarry Smith 
16759bd66eb0SPeter Brune   Input Parameters:
1676ceaaa498SBarry Smith + linesearch  - the linesearch object
16778434afd1SBarry Smith . projectfunc - function for projecting the function to the bounds, see `SNESLineSearchVIProjectFn` for calling sequence
16788434afd1SBarry Smith - normfunc    - function for computing the norm of an active set, see `SNESLineSearchVINormFn` for calling sequence
16799bd66eb0SPeter Brune 
1680f6dfbefdSBarry Smith   Level: advanced
16819bd66eb0SPeter Brune 
1682ceaaa498SBarry Smith   Notes:
168320f4b53cSBarry Smith   The VI solvers require projection of the solution to the feasible set.  `projectfunc` should implement this.
168420f4b53cSBarry Smith 
168520f4b53cSBarry Smith   The VI solvers require special evaluation of the function norm such that the norm is only calculated
168620f4b53cSBarry Smith   on the inactive set.  This should be implemented by `normfunc`.
168720f4b53cSBarry Smith 
16889bcc50f1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetVIFunctions()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()`,
16898434afd1SBarry Smith           `SNESLineSearchVIProjectFn`, `SNESLineSearchVINormFn`
16909bd66eb0SPeter Brune @*/
16918434afd1SBarry Smith PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFn *projectfunc, SNESLineSearchVINormFn *normfunc)
1692d71ae5a4SJacob Faibussowitsch {
16939bd66eb0SPeter Brune   PetscFunctionBegin;
1694f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
16959bd66eb0SPeter Brune   if (projectfunc) linesearch->ops->viproject = projectfunc;
16969bd66eb0SPeter Brune   if (normfunc) linesearch->ops->vinorm = normfunc;
16973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16989bd66eb0SPeter Brune }
16999bd66eb0SPeter Brune 
17009bd66eb0SPeter Brune /*@C
1701f1c6b773SPeter Brune   SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.
17029bd66eb0SPeter Brune 
1703f6dfbefdSBarry Smith   Not Collective
1704f6dfbefdSBarry Smith 
1705f899ff85SJose E. Roman   Input Parameter:
1706f6dfbefdSBarry Smith . linesearch - the line search context, obtain with `SNESGetLineSearch()`
17079bd66eb0SPeter Brune 
17089bd66eb0SPeter Brune   Output Parameters:
17098434afd1SBarry Smith + projectfunc - function for projecting the function to the bounds, see `SNESLineSearchVIProjectFn` for calling sequence
17108434afd1SBarry Smith - normfunc    - function for computing the norm of an active set, see `SNESLineSearchVINormFn ` for calling sequence
17119bd66eb0SPeter Brune 
1712f6dfbefdSBarry Smith   Level: advanced
17139bd66eb0SPeter Brune 
17149bcc50f1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetVIFunctions()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()`,
17158434afd1SBarry Smith           `SNESLineSearchVIProjectFn`, `SNESLineSearchVINormFn`
17169bd66eb0SPeter Brune @*/
17178434afd1SBarry Smith PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFn **projectfunc, SNESLineSearchVINormFn **normfunc)
1718d71ae5a4SJacob Faibussowitsch {
17199bd66eb0SPeter Brune   PetscFunctionBegin;
17209bd66eb0SPeter Brune   if (projectfunc) *projectfunc = linesearch->ops->viproject;
17219bd66eb0SPeter Brune   if (normfunc) *normfunc = linesearch->ops->vinorm;
17223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17239bd66eb0SPeter Brune }
17249bd66eb0SPeter Brune 
1725bf7f4e0aSPeter Brune /*@C
1726420bcc1bSBarry Smith   SNESLineSearchRegister - register a line search type `SNESLineSearchType`
1727ceaaa498SBarry Smith 
1728cc4c1da9SBarry Smith   Logically Collective, No Fortran Support
1729cc4c1da9SBarry Smith 
1730ceaaa498SBarry Smith   Input Parameters:
1731420bcc1bSBarry Smith + sname    - name of the `SNESLineSearchType()`
1732ceaaa498SBarry Smith - function - the creation function for that type
1733ceaaa498SBarry Smith 
1734ceaaa498SBarry Smith   Calling sequence of `function`:
1735ceaaa498SBarry Smith . ls - the line search context
1736bf7f4e0aSPeter Brune 
1737bf7f4e0aSPeter Brune   Level: advanced
1738f6dfbefdSBarry Smith 
1739420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchSetType()`
1740bf7f4e0aSPeter Brune @*/
1741ceaaa498SBarry Smith PetscErrorCode SNESLineSearchRegister(const char sname[], PetscErrorCode (*function)(SNESLineSearch ls))
1742d71ae5a4SJacob Faibussowitsch {
1743bf7f4e0aSPeter Brune   PetscFunctionBegin;
17449566063dSJacob Faibussowitsch   PetscCall(SNESInitializePackage());
17459566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&SNESLineSearchList, sname, function));
17463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1747bf7f4e0aSPeter Brune }
1748