xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision 420bcc1b905230dede3c88f397d1a4e60493adde)
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 
27*420bcc1bSBarry Smith   This does not clear the monitor set with `SNESLineSearchSetDefaultMonitor()` use `SNESLineSearchSetDefaultMonitor`(`ls`,`NULL`) to cancel it
28dcf2fd19SBarry Smith   that one.
29dcf2fd19SBarry Smith 
30*420bcc1bSBarry 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:
56*420bcc1bSBarry Smith   This routine is called by the `SNESLineSearch` implementations.
57dcf2fd19SBarry Smith   It does not typically need to be called by the user.
58dcf2fd19SBarry Smith 
59*420bcc1bSBarry 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
80*420bcc1bSBarry Smith . mctx           - [optional] user-defined context for private data for the monitor routine (use `NULL` if no context is desired)
81*420bcc1bSBarry Smith - monitordestroy - [optional] routine that frees monitor context (may be `NULL`)
82*420bcc1bSBarry Smith 
83*420bcc1bSBarry Smith   Calling sequence of `f`:
84*420bcc1bSBarry Smith + ls   - the `SNESLineSearch` context
85*420bcc1bSBarry Smith - mctx - [optional] user-defined context for private data for the monitor routine
86*420bcc1bSBarry Smith 
87*420bcc1bSBarry Smith   Calling sequence of `monitordestroy`:
88*420bcc1bSBarry Smith . mctx - [optional] user-defined context for private data for the monitor routine
8920f4b53cSBarry Smith 
9020f4b53cSBarry Smith   Level: intermediate
91dcf2fd19SBarry Smith 
92f6dfbefdSBarry Smith   Note:
93dcf2fd19SBarry Smith   Several different monitoring routines may be set by calling
94f6dfbefdSBarry Smith   `SNESLineSearchMonitorSet()` multiple times; all will be called in the
95dcf2fd19SBarry Smith   order in which they were set.
96dcf2fd19SBarry Smith 
97*420bcc1bSBarry Smith   Fortran Note:
98f6dfbefdSBarry Smith   Only a single monitor function can be set for each `SNESLineSearch` object
99dcf2fd19SBarry Smith 
100*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorDefault()`, `SNESLineSearchMonitorCancel()`
101dcf2fd19SBarry Smith @*/
102*420bcc1bSBarry Smith PetscErrorCode SNESLineSearchMonitorSet(SNESLineSearch ls, PetscErrorCode (*f)(SNESLineSearch ls, void *mctx), void *mctx, PetscErrorCode (*monitordestroy)(void **mctx))
103d71ae5a4SJacob Faibussowitsch {
10478064530SBarry Smith   PetscInt  i;
10578064530SBarry Smith   PetscBool identical;
10678064530SBarry Smith 
107dcf2fd19SBarry Smith   PetscFunctionBegin;
108dcf2fd19SBarry Smith   PetscValidHeaderSpecific(ls, SNESLINESEARCH_CLASSID, 1);
10978064530SBarry Smith   for (i = 0; i < ls->numbermonitors; i++) {
1109566063dSJacob Faibussowitsch     PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))f, mctx, monitordestroy, (PetscErrorCode(*)(void))ls->monitorftns[i], ls->monitorcontext[i], ls->monitordestroy[i], &identical));
1113ba16761SJacob Faibussowitsch     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
11278064530SBarry Smith   }
1135f80ce2aSJacob Faibussowitsch   PetscCheck(ls->numbermonitors < MAXSNESLSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many monitors set");
114dcf2fd19SBarry Smith   ls->monitorftns[ls->numbermonitors]      = f;
115dcf2fd19SBarry Smith   ls->monitordestroy[ls->numbermonitors]   = monitordestroy;
116dcf2fd19SBarry Smith   ls->monitorcontext[ls->numbermonitors++] = (void *)mctx;
1173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
118dcf2fd19SBarry Smith }
119dcf2fd19SBarry Smith 
120dcf2fd19SBarry Smith /*@C
121f6dfbefdSBarry Smith   SNESLineSearchMonitorSolutionUpdate - Monitors each update of the function value the linesearch tries
122dcf2fd19SBarry Smith 
123c3339decSBarry Smith   Collective
124dcf2fd19SBarry Smith 
125dcf2fd19SBarry Smith   Input Parameters:
126*420bcc1bSBarry Smith + ls - the `SNESLineSearch` object
127f6dfbefdSBarry Smith - vf - the context for the monitor, in this case it is an `PetscViewerAndFormat`
128dcf2fd19SBarry Smith 
129*420bcc1bSBarry Smith   Options Database Key:
130*420bcc1bSBarry Smith . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine
131dcf2fd19SBarry Smith 
132*420bcc1bSBarry Smith   Level: developer
133*420bcc1bSBarry Smith 
134*420bcc1bSBarry Smith   This is not normally called directly but is passed to `SNESLineSearchMonitorSet()`
135*420bcc1bSBarry Smith 
136*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorSet()`, `SNESMonitorSolution()`
137dcf2fd19SBarry Smith @*/
138d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchMonitorSolutionUpdate(SNESLineSearch ls, PetscViewerAndFormat *vf)
139d71ae5a4SJacob Faibussowitsch {
140d12e167eSBarry Smith   PetscViewer viewer = vf->viewer;
141dcf2fd19SBarry Smith   Vec         Y, W, G;
142dcf2fd19SBarry Smith 
143dcf2fd19SBarry Smith   PetscFunctionBegin;
1449566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetVecs(ls, NULL, NULL, &Y, &W, &G));
1459566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(viewer, vf->format));
1469566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted update to solution \n"));
1479566063dSJacob Faibussowitsch   PetscCall(VecView(Y, viewer));
1489566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted new solution \n"));
1499566063dSJacob Faibussowitsch   PetscCall(VecView(W, viewer));
1509566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted updated function value\n"));
1519566063dSJacob Faibussowitsch   PetscCall(VecView(G, viewer));
1529566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(viewer));
1533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
154dcf2fd19SBarry Smith }
155dcf2fd19SBarry Smith 
156f40b411bSPeter Brune /*@
157*420bcc1bSBarry Smith   SNESLineSearchCreate - Creates a `SNESLineSearch` context.
158f40b411bSPeter Brune 
159f6dfbefdSBarry Smith   Logically Collective
160f40b411bSPeter Brune 
1612fe279fdSBarry Smith   Input Parameter:
162f6dfbefdSBarry Smith . comm - MPI communicator for the line search (typically from the associated `SNES` context).
163f40b411bSPeter Brune 
1642fe279fdSBarry Smith   Output Parameter:
1658e557f58SPeter Brune . outlinesearch - the new line search context
166f40b411bSPeter Brune 
167162e0bf5SPeter Brune   Level: developer
168162e0bf5SPeter Brune 
169f6dfbefdSBarry Smith   Note:
170*420bcc1bSBarry Smith   The preferred calling sequence is to use `SNESGetLineSearch()` to acquire the `SNESLineSearch` instance
171f6dfbefdSBarry Smith   already associated with the `SNES`.
172f40b411bSPeter Brune 
173*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `LineSearchDestroy()`, `SNESGetLineSearch()`
174f40b411bSPeter Brune @*/
175d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch)
176d71ae5a4SJacob Faibussowitsch {
177f1c6b773SPeter Brune   SNESLineSearch linesearch;
178bf388a1fSBarry Smith 
179bf7f4e0aSPeter Brune   PetscFunctionBegin;
1804f572ea9SToby Isaac   PetscAssertPointer(outlinesearch, 2);
1819566063dSJacob Faibussowitsch   PetscCall(SNESInitializePackage());
1820298fd71SBarry Smith   *outlinesearch = NULL;
183f5af7f23SKarl Rupp 
1849566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(linesearch, SNESLINESEARCH_CLASSID, "SNESLineSearch", "Linesearch", "SNESLineSearch", comm, SNESLineSearchDestroy, SNESLineSearchView));
185bf7f4e0aSPeter Brune 
1860298fd71SBarry Smith   linesearch->vec_sol_new  = NULL;
1870298fd71SBarry Smith   linesearch->vec_func_new = NULL;
1880298fd71SBarry Smith   linesearch->vec_sol      = NULL;
1890298fd71SBarry Smith   linesearch->vec_func     = NULL;
1900298fd71SBarry Smith   linesearch->vec_update   = NULL;
1919bd66eb0SPeter Brune 
192bf7f4e0aSPeter Brune   linesearch->lambda       = 1.0;
193bf7f4e0aSPeter Brune   linesearch->fnorm        = 1.0;
194bf7f4e0aSPeter Brune   linesearch->ynorm        = 1.0;
195bf7f4e0aSPeter Brune   linesearch->xnorm        = 1.0;
196422a814eSBarry Smith   linesearch->result       = SNES_LINESEARCH_SUCCEEDED;
197bf7f4e0aSPeter Brune   linesearch->norms        = PETSC_TRUE;
198bf7f4e0aSPeter Brune   linesearch->keeplambda   = PETSC_FALSE;
199bf7f4e0aSPeter Brune   linesearch->damping      = 1.0;
200bf7f4e0aSPeter Brune   linesearch->maxstep      = 1e8;
201bf7f4e0aSPeter Brune   linesearch->steptol      = 1e-12;
202516fe3c3SPeter Brune   linesearch->rtol         = 1e-8;
203516fe3c3SPeter Brune   linesearch->atol         = 1e-15;
204516fe3c3SPeter Brune   linesearch->ltol         = 1e-8;
2050298fd71SBarry Smith   linesearch->precheckctx  = NULL;
2060298fd71SBarry Smith   linesearch->postcheckctx = NULL;
20759405d5eSPeter Brune   linesearch->max_its      = 1;
208bf7f4e0aSPeter Brune   linesearch->setupcalled  = PETSC_FALSE;
2093add74b1SBarry Smith   linesearch->monitor      = NULL;
210bf7f4e0aSPeter Brune   *outlinesearch           = linesearch;
2113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
212bf7f4e0aSPeter Brune }
213bf7f4e0aSPeter Brune 
214f40b411bSPeter Brune /*@
21578bcb3b5SPeter Brune   SNESLineSearchSetUp - Prepares the line search for being applied by allocating
21678bcb3b5SPeter Brune   any required vectors.
217f40b411bSPeter Brune 
218c3339decSBarry Smith   Collective
219f40b411bSPeter Brune 
2202fe279fdSBarry Smith   Input Parameter:
221f6dfbefdSBarry Smith . linesearch - The `SNESLineSearch` instance.
222f40b411bSPeter Brune 
22320f4b53cSBarry Smith   Level: advanced
22420f4b53cSBarry Smith 
225f6dfbefdSBarry Smith   Note:
226f6dfbefdSBarry Smith   For most cases, this needn't be called by users or outside of `SNESLineSearchApply()`.
227cd7522eaSPeter Brune   The only current case where this is called outside of this is for the VI
22878bcb3b5SPeter Brune   solvers, which modify the solution and work vectors before the first call
229f6dfbefdSBarry Smith   of `SNESLineSearchApply()`, requiring the `SNESLineSearch` work vectors to be
230cd7522eaSPeter Brune   allocated upfront.
231cd7522eaSPeter Brune 
232*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchReset()`
233f40b411bSPeter Brune @*/
234d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch)
235d71ae5a4SJacob Faibussowitsch {
236bf7f4e0aSPeter Brune   PetscFunctionBegin;
23748a46eb9SPierre Jolivet   if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC));
238bf7f4e0aSPeter Brune   if (!linesearch->setupcalled) {
23948a46eb9SPierre Jolivet     if (!linesearch->vec_sol_new) PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new));
24048a46eb9SPierre Jolivet     if (!linesearch->vec_func_new) PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new));
241dbbe0bcdSBarry Smith     if (linesearch->ops->setup) PetscUseTypeMethod(linesearch, setup);
2429566063dSJacob Faibussowitsch     if (!linesearch->ops->snesfunc) PetscCall(SNESLineSearchSetFunction(linesearch, SNESComputeFunction));
243bf7f4e0aSPeter Brune     linesearch->lambda      = linesearch->damping;
244bf7f4e0aSPeter Brune     linesearch->setupcalled = PETSC_TRUE;
245bf7f4e0aSPeter Brune   }
2463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
247bf7f4e0aSPeter Brune }
248bf7f4e0aSPeter Brune 
249f40b411bSPeter Brune /*@
250f6dfbefdSBarry Smith   SNESLineSearchReset - Undoes the `SNESLineSearchSetUp()` and deletes any `Vec`s or `Mat`s allocated by the line search.
251f40b411bSPeter Brune 
252c3339decSBarry Smith   Collective
253f40b411bSPeter Brune 
2542fe279fdSBarry Smith   Input Parameter:
255f6dfbefdSBarry Smith . linesearch - The `SNESLineSearch` instance.
256f40b411bSPeter Brune 
25720f4b53cSBarry Smith   Level: developer
25820f4b53cSBarry Smith 
259f6dfbefdSBarry Smith   Note:
260f6dfbefdSBarry Smith   Usually only called by `SNESReset()`
261f190f2fcSBarry Smith 
262*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchSetUp()`
263f40b411bSPeter Brune @*/
264d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch)
265d71ae5a4SJacob Faibussowitsch {
266bf7f4e0aSPeter Brune   PetscFunctionBegin;
267dbbe0bcdSBarry Smith   if (linesearch->ops->reset) PetscUseTypeMethod(linesearch, reset);
268f5af7f23SKarl Rupp 
2699566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&linesearch->vec_sol_new));
2709566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&linesearch->vec_func_new));
271bf7f4e0aSPeter Brune 
2729566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(linesearch->nwork, &linesearch->work));
273f5af7f23SKarl Rupp 
274bf7f4e0aSPeter Brune   linesearch->nwork       = 0;
275bf7f4e0aSPeter Brune   linesearch->setupcalled = PETSC_FALSE;
2763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
277bf7f4e0aSPeter Brune }
278bf7f4e0aSPeter Brune 
279ed07d7d7SPeter Brune /*@C
280f6dfbefdSBarry Smith   SNESLineSearchSetFunction - Sets the function evaluation used by the `SNES` line search
281f6dfbefdSBarry Smith   `
282e4094ef1SJacob Faibussowitsch 
283ed07d7d7SPeter Brune   Input Parameters:
284e4094ef1SJacob Faibussowitsch + linesearch - the `SNESLineSearch` context
285e4094ef1SJacob Faibussowitsch - func       - function evaluation routine, this is usually the function provided with `SNESSetFunction()`
286ed07d7d7SPeter Brune 
287*420bcc1bSBarry Smith   Calling sequence of `func`:
288*420bcc1bSBarry Smith + snes - the `SNES` with which the `SNESLineSearch` context is associated with
289*420bcc1bSBarry Smith . x    - the input vector
290*420bcc1bSBarry Smith - f    - the computed value of the function
291*420bcc1bSBarry Smith 
292ed07d7d7SPeter Brune   Level: developer
293ed07d7d7SPeter Brune 
294*420bcc1bSBarry Smith   Note:
295*420bcc1bSBarry Smith   By default the `SNESLineSearch` uses the function provided by `SNESSetFunction()` so this is rarely needed
296*420bcc1bSBarry Smith 
297*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESSetFunction()`
298ed07d7d7SPeter Brune @*/
299*420bcc1bSBarry Smith PetscErrorCode SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES snes, Vec x, Vec f))
300d71ae5a4SJacob Faibussowitsch {
301ed07d7d7SPeter Brune   PetscFunctionBegin;
302ed07d7d7SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
303ed07d7d7SPeter Brune   linesearch->ops->snesfunc = func;
3043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
305ed07d7d7SPeter Brune }
306ed07d7d7SPeter Brune 
30786d74e61SPeter Brune /*@C
308*420bcc1bSBarry Smith   SNESLineSearchSetPreCheck - Sets a function that is called after the initial search direction has been computed but
309*420bcc1bSBarry Smith   before the line search routine has been applied. Allows adjusting the result of (usually a linear solve) that
310f190f2fcSBarry Smith   determined the search direction.
31186d74e61SPeter Brune 
312c3339decSBarry Smith   Logically Collective
31386d74e61SPeter Brune 
31486d74e61SPeter Brune   Input Parameters:
315f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context
316*420bcc1bSBarry Smith . func       - [optional] function evaluation routine
31720f4b53cSBarry Smith - ctx        - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
31886d74e61SPeter Brune 
319*420bcc1bSBarry Smith   Calling sequence of `func`:
320*420bcc1bSBarry Smith + ls        - the `SNESLineSearch` context
321*420bcc1bSBarry Smith . x         - the current solution
322*420bcc1bSBarry Smith . d         - the current seach direction
323*420bcc1bSBarry Smith . changed_d - indicates if the search direction has been changed
324*420bcc1bSBarry Smith - ctx       - the context passed to `SNESLineSearchSetPreCheck()`
325*420bcc1bSBarry Smith 
32686d74e61SPeter Brune   Level: intermediate
32786d74e61SPeter Brune 
328f6dfbefdSBarry Smith   Note:
329*420bcc1bSBarry Smith   Use `SNESLineSearchSetPostCheck()` to change the step after the line search is complete.
330f0b84518SBarry Smith 
331f0b84518SBarry Smith   Use `SNESVISetVariableBounds()` and `SNESVISetComputeVariableBounds()` to cause `SNES` to automatically control the ranges of variables allowed.
332f0b84518SBarry Smith 
333*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()`,
334869ba2dcSBarry Smith           `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`
335f0b84518SBarry Smith 
33686d74e61SPeter Brune @*/
337*420bcc1bSBarry Smith PetscErrorCode SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch ls, Vec x, Vec d, PetscBool *changed_d, void *ctx), void *ctx)
338d71ae5a4SJacob Faibussowitsch {
3399bd66eb0SPeter Brune   PetscFunctionBegin;
340f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
341f190f2fcSBarry Smith   if (func) linesearch->ops->precheck = func;
34286d74e61SPeter Brune   if (ctx) linesearch->precheckctx = ctx;
3433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
34486d74e61SPeter Brune }
34586d74e61SPeter Brune 
34686d74e61SPeter Brune /*@C
34753deb899SBarry Smith   SNESLineSearchGetPreCheck - Gets the pre-check function for the line search routine.
34886d74e61SPeter Brune 
349f899ff85SJose E. Roman   Input Parameter:
350f6dfbefdSBarry Smith . linesearch - the `SNESLineSearch` context
35186d74e61SPeter Brune 
35286d74e61SPeter Brune   Output Parameters:
353*420bcc1bSBarry Smith + func - [optional] function evaluation routine,  for calling sequence see `SNESLineSearchSetPreCheck()`
35420f4b53cSBarry Smith - ctx  - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
35586d74e61SPeter Brune 
35686d74e61SPeter Brune   Level: intermediate
35786d74e61SPeter Brune 
358*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchSetPostCheck()`
35986d74e61SPeter Brune @*/
360d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch, Vec, Vec, PetscBool *, void *), void **ctx)
361d71ae5a4SJacob Faibussowitsch {
3629bd66eb0SPeter Brune   PetscFunctionBegin;
363f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
364f190f2fcSBarry Smith   if (func) *func = linesearch->ops->precheck;
36586d74e61SPeter Brune   if (ctx) *ctx = linesearch->precheckctx;
3663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
36786d74e61SPeter Brune }
36886d74e61SPeter Brune 
36986d74e61SPeter Brune /*@C
370f190f2fcSBarry Smith   SNESLineSearchSetPostCheck - Sets a user function that is called after the line search has been applied to determine the step
371f190f2fcSBarry Smith   direction and length. Allows the user a chance to change or override the decision of the line search routine
37286d74e61SPeter Brune 
37320f4b53cSBarry Smith   Logically Collective
37486d74e61SPeter Brune 
37586d74e61SPeter Brune   Input Parameters:
376f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context
377*420bcc1bSBarry Smith . func       - [optional] function evaluation routine
37820f4b53cSBarry Smith - ctx        - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
37986d74e61SPeter Brune 
380*420bcc1bSBarry Smith   Calling sequence of `func`:
381*420bcc1bSBarry Smith + ls        - the `SNESLineSearch` context
382*420bcc1bSBarry Smith . x         - the current solution
383*420bcc1bSBarry Smith . d         - the current seach direction
384*420bcc1bSBarry Smith . w         - $ w = x + lamba*d $ for some lambda
385*420bcc1bSBarry Smith . changed_d - indicates if the search direction `d` has been changed
386*420bcc1bSBarry Smith . changed_w - indicates `w` has been changed
387*420bcc1bSBarry Smith - ctx       - the context passed to `SNESLineSearchSetPreCheck()`
388*420bcc1bSBarry Smith 
38986d74e61SPeter Brune   Level: intermediate
39086d74e61SPeter Brune 
391f0b84518SBarry Smith   Notes:
392*420bcc1bSBarry Smith   Use `SNESLineSearchSetPreCheck()` to change the step before the line search is complete.
393f0b84518SBarry Smith 
394f0b84518SBarry Smith   Use `SNESVISetVariableBounds()` and `SNESVISetComputeVariableBounds()` to cause `SNES` to automatically control the ranges of variables allowed.
395f0b84518SBarry Smith 
396*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchGetPostCheck()`,
397e4094ef1SJacob Faibussowitsch           `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`
39886d74e61SPeter Brune @*/
399*420bcc1bSBarry 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)
400d71ae5a4SJacob Faibussowitsch {
40186d74e61SPeter Brune   PetscFunctionBegin;
402f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
403f190f2fcSBarry Smith   if (func) linesearch->ops->postcheck = func;
40486d74e61SPeter Brune   if (ctx) linesearch->postcheckctx = ctx;
4053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
40686d74e61SPeter Brune }
40786d74e61SPeter Brune 
40886d74e61SPeter Brune /*@C
409f1c6b773SPeter Brune   SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine.
41086d74e61SPeter Brune 
411f899ff85SJose E. Roman   Input Parameter:
412f6dfbefdSBarry Smith . linesearch - the `SNESLineSearch` context
41386d74e61SPeter Brune 
41486d74e61SPeter Brune   Output Parameters:
415*420bcc1bSBarry Smith + func - [optional] function evaluation routine, see for the calling sequence `SNESLineSearchSetPostCheck()`
41620f4b53cSBarry Smith - ctx  - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
41786d74e61SPeter Brune 
41886d74e61SPeter Brune   Level: intermediate
41986d74e61SPeter Brune 
420*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`
42186d74e61SPeter Brune @*/
422d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx)
423d71ae5a4SJacob Faibussowitsch {
4249bd66eb0SPeter Brune   PetscFunctionBegin;
425f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
426f190f2fcSBarry Smith   if (func) *func = linesearch->ops->postcheck;
42786d74e61SPeter Brune   if (ctx) *ctx = linesearch->postcheckctx;
4283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
42986d74e61SPeter Brune }
43086d74e61SPeter Brune 
431f40b411bSPeter Brune /*@
432f1c6b773SPeter Brune   SNESLineSearchPreCheck - Prepares the line search for being applied.
433f40b411bSPeter Brune 
434c3339decSBarry Smith   Logically Collective
435f40b411bSPeter Brune 
436f40b411bSPeter Brune   Input Parameters:
4377b1df9c1SPeter Brune + linesearch - The linesearch instance.
4387b1df9c1SPeter Brune . X          - The current solution
4397b1df9c1SPeter Brune - Y          - The step direction
440f40b411bSPeter Brune 
4412fe279fdSBarry Smith   Output Parameter:
442*420bcc1bSBarry Smith . changed - Indicator that the precheck routine has changed `Y`
443f40b411bSPeter Brune 
444f0b84518SBarry Smith   Level: advanced
445f40b411bSPeter Brune 
446f6dfbefdSBarry Smith   Note:
447*420bcc1bSBarry Smith   This calls any function provided with `SNESLineSearchSetPreCheck()` and is called automatically inside the line search routines
448f6dfbefdSBarry Smith 
449*420bcc1bSBarry Smith   Developer Note:
450*420bcc1bSBarry Smith   The use of `PetscObjectGetState()` would eliminate the need for the `changed` argument to be provided
451*420bcc1bSBarry Smith 
452*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`,
453fe8e7dddSPierre Jolivet           `SNESLineSearchGetPostCheck()`
454f40b411bSPeter Brune @*/
455d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed)
456d71ae5a4SJacob Faibussowitsch {
457bf7f4e0aSPeter Brune   PetscFunctionBegin;
458bf7f4e0aSPeter Brune   *changed = PETSC_FALSE;
4596b2b7091SBarry Smith   if (linesearch->ops->precheck) {
460dbbe0bcdSBarry Smith     PetscUseTypeMethod(linesearch, precheck, X, Y, changed, linesearch->precheckctx);
46138bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch, *changed, 4);
462bf7f4e0aSPeter Brune   }
4633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
464bf7f4e0aSPeter Brune }
465bf7f4e0aSPeter Brune 
466f40b411bSPeter Brune /*@
467ef46b1a6SStefano Zampini   SNESLineSearchPostCheck - Hook to modify step direction or updated solution after a successful linesearch
468f40b411bSPeter Brune 
469c3339decSBarry Smith   Logically Collective
470f40b411bSPeter Brune 
471f40b411bSPeter Brune   Input Parameters:
4727b1df9c1SPeter Brune + linesearch - The line search context
4737b1df9c1SPeter Brune . X          - The last solution
4747b1df9c1SPeter Brune . Y          - The step direction
475*420bcc1bSBarry Smith - W          - The updated solution, $ W = X + lambda*Y $ for some lambda
476f40b411bSPeter Brune 
477f40b411bSPeter Brune   Output Parameters:
478*420bcc1bSBarry Smith + changed_Y - Indicator if the direction `Y` has been changed.
479*420bcc1bSBarry Smith - changed_W - Indicator if the new candidate solution `W` has been changed.
480f40b411bSPeter Brune 
48120f4b53cSBarry Smith   Level: developer
48220f4b53cSBarry Smith 
483f6dfbefdSBarry Smith   Note:
484*420bcc1bSBarry Smith   This calls any function provided with `SNESLineSearchSetPostCheck()` and is called automatically inside the line search routines
485f6dfbefdSBarry Smith 
486*420bcc1bSBarry Smith   Developer Note:
487*420bcc1bSBarry Smith   The use of `PetscObjectGetState()` would eliminate the need for the `changed_Y` and `changed_W` arguments to be provided
488*420bcc1bSBarry Smith 
489*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPrecheck()`, `SNESLineSearchGetPrecheck()`
490f40b411bSPeter Brune @*/
491d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W)
492d71ae5a4SJacob Faibussowitsch {
493bf7f4e0aSPeter Brune   PetscFunctionBegin;
494bf7f4e0aSPeter Brune   *changed_Y = PETSC_FALSE;
495bf7f4e0aSPeter Brune   *changed_W = PETSC_FALSE;
4966b2b7091SBarry Smith   if (linesearch->ops->postcheck) {
497dbbe0bcdSBarry Smith     PetscUseTypeMethod(linesearch, postcheck, X, Y, W, changed_Y, changed_W, linesearch->postcheckctx);
49838bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch, *changed_Y, 5);
49938bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch, *changed_W, 6);
50086d74e61SPeter Brune   }
5013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
50286d74e61SPeter Brune }
50386d74e61SPeter Brune 
50486d74e61SPeter Brune /*@C
50586d74e61SPeter Brune   SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration
50686d74e61SPeter Brune 
507c3339decSBarry Smith   Logically Collective
50886d74e61SPeter Brune 
5094165533cSJose E. Roman   Input Parameters:
510*420bcc1bSBarry Smith + linesearch - the line search context
51186d74e61SPeter Brune . X          - base state for this step
512907376e6SBarry Smith - ctx        - context for this function
51386d74e61SPeter Brune 
51497bb3fdcSJose E. Roman   Input/Output Parameter:
51597bb3fdcSJose E. Roman . Y - correction, possibly modified
51697bb3fdcSJose E. Roman 
51797bb3fdcSJose E. Roman   Output Parameter:
518*420bcc1bSBarry Smith . changed - flag indicating that `Y` was modified
51986d74e61SPeter Brune 
520*420bcc1bSBarry Smith   Options Database Keys:
521cd7522eaSPeter Brune + -snes_linesearch_precheck_picard       - activate this routine
522cd7522eaSPeter Brune - -snes_linesearch_precheck_picard_angle - angle
52386d74e61SPeter Brune 
52486d74e61SPeter Brune   Level: advanced
52586d74e61SPeter Brune 
52686d74e61SPeter Brune   Notes:
527f6dfbefdSBarry Smith   This function should be passed to `SNESLineSearchSetPreCheck()`
52886d74e61SPeter Brune 
52986d74e61SPeter Brune   The justification for this method involves the linear convergence of a Picard iteration
53086d74e61SPeter Brune   so the Picard linearization should be provided in place of the "Jacobian". This correction
53186d74e61SPeter Brune   is generally not useful when using a Newton linearization.
53286d74e61SPeter Brune 
533e4094ef1SJacob Faibussowitsch   References:
534f6dfbefdSBarry Smith   . - * - Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology.
53586d74e61SPeter Brune 
536*420bcc1bSBarry Smith   Developer Note:
537*420bcc1bSBarry Smith   The use of `PetscObjectGetState()` would eliminate the need for the `changed` argument to be provided
538*420bcc1bSBarry Smith 
539*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESSetPicard()`, `SNESGetLineSearch()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchSetPostCheck()`
54086d74e61SPeter Brune @*/
541d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed, void *ctx)
542d71ae5a4SJacob Faibussowitsch {
54386d74e61SPeter Brune   PetscReal   angle = *(PetscReal *)linesearch->precheckctx;
54486d74e61SPeter Brune   Vec         Ylast;
54586d74e61SPeter Brune   PetscScalar dot;
54686d74e61SPeter Brune   PetscInt    iter;
54786d74e61SPeter Brune   PetscReal   ynorm, ylastnorm, theta, angle_radians;
54886d74e61SPeter Brune   SNES        snes;
54986d74e61SPeter Brune 
55086d74e61SPeter Brune   PetscFunctionBegin;
5519566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
5529566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)snes, "SNESLineSearchPreCheckPicard_Ylast", (PetscObject *)&Ylast));
55386d74e61SPeter Brune   if (!Ylast) {
5549566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(Y, &Ylast));
5559566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)snes, "SNESLineSearchPreCheckPicard_Ylast", (PetscObject)Ylast));
5569566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)Ylast));
55786d74e61SPeter Brune   }
5589566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &iter));
55986d74e61SPeter Brune   if (iter < 2) {
5609566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y, Ylast));
56186d74e61SPeter Brune     *changed = PETSC_FALSE;
5623ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
56386d74e61SPeter Brune   }
56486d74e61SPeter Brune 
5659566063dSJacob Faibussowitsch   PetscCall(VecDot(Y, Ylast, &dot));
5669566063dSJacob Faibussowitsch   PetscCall(VecNorm(Y, NORM_2, &ynorm));
5679566063dSJacob Faibussowitsch   PetscCall(VecNorm(Ylast, NORM_2, &ylastnorm));
56886d74e61SPeter Brune   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
569255453a1SBarry Smith   theta         = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm), -1.0, 1.0));
57086d74e61SPeter Brune   angle_radians = angle * PETSC_PI / 180.;
57186d74e61SPeter Brune   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
57286d74e61SPeter Brune     /* Modify the step Y */
57386d74e61SPeter Brune     PetscReal alpha, ydiffnorm;
5749566063dSJacob Faibussowitsch     PetscCall(VecAXPY(Ylast, -1.0, Y));
5759566063dSJacob Faibussowitsch     PetscCall(VecNorm(Ylast, NORM_2, &ydiffnorm));
576f85e2ce2SBarry Smith     alpha = (ydiffnorm > .001 * ylastnorm) ? ylastnorm / ydiffnorm : 1000.0;
5779566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y, Ylast));
5789566063dSJacob Faibussowitsch     PetscCall(VecScale(Y, alpha));
5799566063dSJacob 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));
580a47ec85fSBarry Smith     *changed = PETSC_TRUE;
58186d74e61SPeter Brune   } else {
5829566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Angle %14.12e degrees exceeds threshold %14.12e, no correction applied\n", (double)(theta * 180. / PETSC_PI), (double)angle));
5839566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y, Ylast));
58486d74e61SPeter Brune     *changed = PETSC_FALSE;
585bf7f4e0aSPeter Brune   }
5863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
587bf7f4e0aSPeter Brune }
588bf7f4e0aSPeter Brune 
589f40b411bSPeter Brune /*@
590cd7522eaSPeter Brune   SNESLineSearchApply - Computes the line-search update.
591f40b411bSPeter Brune 
592c3339decSBarry Smith   Collective
593f40b411bSPeter Brune 
594f40b411bSPeter Brune   Input Parameters:
5958e557f58SPeter Brune + linesearch - The line search context
5968e557f58SPeter Brune - Y          - The search direction
597f40b411bSPeter Brune 
5986b867d5aSJose E. Roman   Input/Output Parameters:
5996b867d5aSJose E. Roman + X     - The current solution, on output the new solution
600*420bcc1bSBarry Smith . F     - The current function value, on output the new function value at the solution value `X`
601*420bcc1bSBarry Smith - fnorm - The current norm of `F`, on output the new norm of `F`
602f40b411bSPeter Brune 
603cd7522eaSPeter Brune   Options Database Keys:
6040b00b554SBarry Smith + -snes_linesearch_type                - basic (or equivalently none), bt, l2, cp, nleqerr, shell
605dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches
6061fe24845SBarry Smith . -snes_linesearch_damping             - The linesearch damping parameter, default is 1.0 (no damping)
6071fe24845SBarry Smith . -snes_linesearch_norms               - Turn on/off the linesearch norms computation (SNESLineSearchSetComputeNorms())
6083c7d6663SPeter Brune . -snes_linesearch_keeplambda          - Keep the previous search length as the initial guess
6093c7d6663SPeter Brune - -snes_linesearch_max_it              - The number of iterations for iterative line searches
610cd7522eaSPeter Brune 
611e4094ef1SJacob Faibussowitsch   Level: intermediate
61220f4b53cSBarry Smith 
613cd7522eaSPeter Brune   Notes:
614f6dfbefdSBarry Smith   This is typically called from within a `SNESSolve()` implementation in order to
615f6dfbefdSBarry Smith   help with convergence of the nonlinear method.  Various `SNES` types use line searches
616cd7522eaSPeter Brune   in different ways, but the overarching theme is that a line search is used to determine
617cd7522eaSPeter Brune   an optimal damping parameter of a step at each iteration of the method.  Each
618f6dfbefdSBarry Smith   application of the line search may invoke `SNESComputeFunction()` several times, and
619cd7522eaSPeter Brune   therefore may be fairly expensive.
620cd7522eaSPeter Brune 
621*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchPreCheck()`, `SNESLineSearchPostCheck()`, `SNESSolve()`, `SNESComputeFunction()`, `SNESLineSearchSetComputeNorms()`,
622db781477SPatrick Sanan           `SNESLineSearchType`, `SNESLineSearchSetType()`
623f40b411bSPeter Brune @*/
624d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal *fnorm, Vec Y)
625d71ae5a4SJacob Faibussowitsch {
626bf388a1fSBarry Smith   PetscFunctionBegin;
627f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
628bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
629bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
630064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(Y, VEC_CLASSID, 5);
631bf7f4e0aSPeter Brune 
632422a814eSBarry Smith   linesearch->result = SNES_LINESEARCH_SUCCEEDED;
633bf7f4e0aSPeter Brune 
634bf7f4e0aSPeter Brune   linesearch->vec_sol    = X;
635bf7f4e0aSPeter Brune   linesearch->vec_update = Y;
636bf7f4e0aSPeter Brune   linesearch->vec_func   = F;
637bf7f4e0aSPeter Brune 
6389566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchSetUp(linesearch));
639bf7f4e0aSPeter Brune 
640f5af7f23SKarl Rupp   if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */
641bf7f4e0aSPeter Brune 
642f5af7f23SKarl Rupp   if (fnorm) linesearch->fnorm = *fnorm;
64348a46eb9SPierre Jolivet   else PetscCall(VecNorm(F, NORM_2, &linesearch->fnorm));
644bf7f4e0aSPeter Brune 
6459566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(SNESLINESEARCH_Apply, linesearch, X, F, Y));
646bf7f4e0aSPeter Brune 
647dbbe0bcdSBarry Smith   PetscUseTypeMethod(linesearch, apply);
648bf7f4e0aSPeter Brune 
6499566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(SNESLINESEARCH_Apply, linesearch, X, F, Y));
650bf7f4e0aSPeter Brune 
651f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
6523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
653bf7f4e0aSPeter Brune }
654bf7f4e0aSPeter Brune 
655f40b411bSPeter Brune /*@
656f1c6b773SPeter Brune   SNESLineSearchDestroy - Destroys the line search instance.
657f40b411bSPeter Brune 
658c3339decSBarry Smith   Collective
659f40b411bSPeter Brune 
6602fe279fdSBarry Smith   Input Parameter:
6618e557f58SPeter Brune . linesearch - The line search context
662f40b411bSPeter Brune 
66384238204SBarry Smith   Level: developer
664f40b411bSPeter Brune 
665*420bcc1bSBarry Smith   Note:
666*420bcc1bSBarry Smith   The line search in `SNES` is automatically called on `SNESDestroy()` so this call is rarely needed
667*420bcc1bSBarry Smith 
668*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchReset()`, `SNESDestroy()`
669f40b411bSPeter Brune @*/
670d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchDestroy(SNESLineSearch *linesearch)
671d71ae5a4SJacob Faibussowitsch {
672bf7f4e0aSPeter Brune   PetscFunctionBegin;
6733ba16761SJacob Faibussowitsch   if (!*linesearch) PetscFunctionReturn(PETSC_SUCCESS);
674f1c6b773SPeter Brune   PetscValidHeaderSpecific((*linesearch), SNESLINESEARCH_CLASSID, 1);
6759371c9d4SSatish Balay   if (--((PetscObject)(*linesearch))->refct > 0) {
6769371c9d4SSatish Balay     *linesearch = NULL;
6773ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
6789371c9d4SSatish Balay   }
6799566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsViewOff((PetscObject)*linesearch));
6809566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchReset(*linesearch));
6813ba16761SJacob Faibussowitsch   PetscTryTypeMethod(*linesearch, destroy);
6829566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&(*linesearch)->monitor));
6839566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchMonitorCancel((*linesearch)));
6849566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(linesearch));
6853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
686bf7f4e0aSPeter Brune }
687bf7f4e0aSPeter Brune 
688f40b411bSPeter Brune /*@
689dcf2fd19SBarry Smith   SNESLineSearchSetDefaultMonitor - Turns on/off printing useful information and debugging output about the line search.
690bf7f4e0aSPeter Brune 
691c3339decSBarry Smith   Logically Collective
692f6dfbefdSBarry Smith 
693bf7f4e0aSPeter Brune   Input Parameters:
694dcf2fd19SBarry Smith + linesearch - the linesearch object
69520f4b53cSBarry Smith - viewer     - an `PETSCVIEWERASCII` `PetscViewer` or `NULL` to turn off monitor
696bf7f4e0aSPeter Brune 
697f6dfbefdSBarry Smith   Options Database Key:
698dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - enables the monitor
699bf7f4e0aSPeter Brune 
700bf7f4e0aSPeter Brune   Level: intermediate
701bf7f4e0aSPeter Brune 
702e4094ef1SJacob Faibussowitsch   Developer Notes:
703f6dfbefdSBarry Smith   This monitor is implemented differently than the other line search monitors that are set with
704f6dfbefdSBarry Smith   `SNESLineSearchMonitorSet()` since it is called in many locations of the line search routines to display aspects of the
705d12e167eSBarry Smith   line search that are not visible to the other monitors.
706dcf2fd19SBarry Smith 
707*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `PETSCVIEWERASCII`, `SNESGetLineSearch()`, `SNESLineSearchGetDefaultMonitor()`, `PetscViewer`, `SNESLineSearchSetMonitor()`,
708f6dfbefdSBarry Smith           `SNESLineSearchMonitorSetFromOptions()`
709bf7f4e0aSPeter Brune @*/
710d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetDefaultMonitor(SNESLineSearch linesearch, PetscViewer viewer)
711d71ae5a4SJacob Faibussowitsch {
712bf7f4e0aSPeter Brune   PetscFunctionBegin;
7139566063dSJacob Faibussowitsch   if (viewer) PetscCall(PetscObjectReference((PetscObject)viewer));
7149566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&linesearch->monitor));
715dcf2fd19SBarry Smith   linesearch->monitor = viewer;
7163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
717bf7f4e0aSPeter Brune }
718bf7f4e0aSPeter Brune 
719f40b411bSPeter Brune /*@
720*420bcc1bSBarry Smith   SNESLineSearchGetDefaultMonitor - Gets the `PetscViewer` instance for the default line search monitor that is turned on with `SNESLineSearchSetDefaultMonitor()`
721f6dfbefdSBarry Smith 
722c3339decSBarry Smith   Logically Collective
7236a388c36SPeter Brune 
724f190f2fcSBarry Smith   Input Parameter:
725*420bcc1bSBarry Smith . linesearch - the line search context
726f40b411bSPeter Brune 
727f190f2fcSBarry Smith   Output Parameter:
7288e557f58SPeter Brune . monitor - monitor context
729f40b411bSPeter Brune 
730f40b411bSPeter Brune   Level: intermediate
731f40b411bSPeter Brune 
732*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchSetDefaultMonitor()`, `PetscViewer`
733f40b411bSPeter Brune @*/
734d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetDefaultMonitor(SNESLineSearch linesearch, PetscViewer *monitor)
735d71ae5a4SJacob Faibussowitsch {
7366a388c36SPeter Brune   PetscFunctionBegin;
737f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
7386a388c36SPeter Brune   *monitor = linesearch->monitor;
7393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7406a388c36SPeter Brune }
7416a388c36SPeter Brune 
742dcf2fd19SBarry Smith /*@C
743*420bcc1bSBarry Smith   SNESLineSearchMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated in the options database
744dcf2fd19SBarry Smith 
745c3339decSBarry Smith   Collective
746dcf2fd19SBarry Smith 
747dcf2fd19SBarry Smith   Input Parameters:
748*420bcc1bSBarry Smith + ls           - `SNESLineSearch` object to monitor
749*420bcc1bSBarry Smith . name         - the monitor type
750dcf2fd19SBarry Smith . help         - message indicating what monitoring is done
751dcf2fd19SBarry Smith . manual       - manual page for the monitor
752dcf2fd19SBarry Smith . monitor      - the monitor function
753f6dfbefdSBarry 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`
754dcf2fd19SBarry Smith 
755*420bcc1bSBarry Smith   Calling sequence of `monitor`:
756*420bcc1bSBarry Smith + ls - `SNESLineSearch` object being monitored
757*420bcc1bSBarry Smith - vf - a `PetscViewerAndFormat` struct that provides the `PetscViewer` and `PetscViewerFormat` being used
758dcf2fd19SBarry Smith 
759*420bcc1bSBarry Smith   Calling sequence of `monitorsetup`:
760*420bcc1bSBarry Smith + ls - `SNESLineSearch` object being monitored
761*420bcc1bSBarry Smith - vf - a `PetscViewerAndFormat` struct that provides the `PetscViewer` and `PetscViewerFormat` being used
762*420bcc1bSBarry Smith 
763*420bcc1bSBarry Smith   Level: advanced
764*420bcc1bSBarry Smith 
765*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetMonitor()`, `PetscOptionsGetViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
766db781477SPatrick Sanan           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
767e4094ef1SJacob Faibussowitsch           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
768db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
769c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
770db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
771db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
772dcf2fd19SBarry Smith @*/
773*420bcc1bSBarry 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))
774d71ae5a4SJacob Faibussowitsch {
775dcf2fd19SBarry Smith   PetscViewer       viewer;
776dcf2fd19SBarry Smith   PetscViewerFormat format;
777dcf2fd19SBarry Smith   PetscBool         flg;
778dcf2fd19SBarry Smith 
779dcf2fd19SBarry Smith   PetscFunctionBegin;
7809566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)ls), ((PetscObject)ls)->options, ((PetscObject)ls)->prefix, name, &viewer, &format, &flg));
781dcf2fd19SBarry Smith   if (flg) {
782d12e167eSBarry Smith     PetscViewerAndFormat *vf;
7839566063dSJacob Faibussowitsch     PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf));
7849566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)viewer));
7851baa6e33SBarry Smith     if (monitorsetup) PetscCall((*monitorsetup)(ls, vf));
7869566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchMonitorSet(ls, (PetscErrorCode(*)(SNESLineSearch, void *))monitor, vf, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy));
787dcf2fd19SBarry Smith   }
7883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
789dcf2fd19SBarry Smith }
790dcf2fd19SBarry Smith 
791f40b411bSPeter Brune /*@
792f1c6b773SPeter Brune   SNESLineSearchSetFromOptions - Sets options for the line search
793f40b411bSPeter Brune 
794c3339decSBarry Smith   Logically Collective
795f6dfbefdSBarry Smith 
796f899ff85SJose E. Roman   Input Parameter:
797*420bcc1bSBarry Smith . linesearch - a `SNESLineSearch` line search context
798f40b411bSPeter Brune 
799f40b411bSPeter Brune   Options Database Keys:
8000b00b554SBarry Smith + -snes_linesearch_type <type>                                      - basic (or equivalently none), bt, l2, cp, nleqerr, shell
8013c7d6663SPeter Brune . -snes_linesearch_order <order>                                    - 1, 2, 3.  Most types only support certain orders (bt supports 2 or 3)
802f6dfbefdSBarry Smith . -snes_linesearch_norms                                            - Turn on/off the linesearch norms for the basic linesearch typem (`SNESLineSearchSetComputeNorms()`)
80371eef1aeSPeter Brune . -snes_linesearch_minlambda                                        - The minimum step length
8041a9b3a06SPeter Brune . -snes_linesearch_maxstep                                          - The maximum step size
8051a9b3a06SPeter Brune . -snes_linesearch_rtol                                             - Relative tolerance for iterative line searches
8061a9b3a06SPeter Brune . -snes_linesearch_atol                                             - Absolute tolerance for iterative line searches
8071a9b3a06SPeter Brune . -snes_linesearch_ltol                                             - Change in lambda tolerance for iterative line searches
8081a9b3a06SPeter Brune . -snes_linesearch_max_it                                           - The number of iterations for iterative line searches
809dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename]                              - Print progress of line searches
810dcf2fd19SBarry Smith . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine
8118e557f58SPeter Brune . -snes_linesearch_damping                                          - The linesearch damping parameter
812cd7522eaSPeter Brune . -snes_linesearch_keeplambda                                       - Keep the previous search length as the initial guess.
8131a9b3a06SPeter Brune . -snes_linesearch_precheck_picard                                  - Use precheck that speeds up convergence of picard method
814d8d34be6SBarry Smith - -snes_linesearch_precheck_picard_angle                            - Angle used in Picard precheck method
815f40b411bSPeter Brune 
816f40b411bSPeter Brune   Level: intermediate
817f40b411bSPeter Brune 
818*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchCreate()`, `SNESLineSearchSetOrder()`, `SNESLineSearchSetType()`, `SNESLineSearchSetTolerances()`, `SNESLineSearchSetDamping()`, `SNESLineSearchPreCheckPicard()`,
819db781477SPatrick Sanan           `SNESLineSearchType`, `SNESLineSearchSetComputeNorms()`
820f40b411bSPeter Brune @*/
821d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch)
822d71ae5a4SJacob Faibussowitsch {
8231a4f838cSPeter Brune   const char *deft = SNESLINESEARCHBASIC;
824bf7f4e0aSPeter Brune   char        type[256];
825bf7f4e0aSPeter Brune   PetscBool   flg, set;
826dcf2fd19SBarry Smith   PetscViewer viewer;
827bf388a1fSBarry Smith 
828bf7f4e0aSPeter Brune   PetscFunctionBegin;
8299566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchRegisterAll());
830bf7f4e0aSPeter Brune 
831d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)linesearch);
832f5af7f23SKarl Rupp   if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name;
8339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-snes_linesearch_type", "Linesearch type", "SNESLineSearchSetType", SNESLineSearchList, deft, type, 256, &flg));
834bf7f4e0aSPeter Brune   if (flg) {
8359566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetType(linesearch, type));
836bf7f4e0aSPeter Brune   } else if (!((PetscObject)linesearch)->type_name) {
8379566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetType(linesearch, deft));
838bf7f4e0aSPeter Brune   }
839bf7f4e0aSPeter Brune 
8409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)linesearch), ((PetscObject)linesearch)->options, ((PetscObject)linesearch)->prefix, "-snes_linesearch_monitor", &viewer, NULL, &set));
841dcf2fd19SBarry Smith   if (set) {
8429566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetDefaultMonitor(linesearch, viewer));
8439566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
844dcf2fd19SBarry Smith   }
8459566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchMonitorSetFromOptions(linesearch, "-snes_linesearch_monitor_solution_update", "View correction at each iteration", "SNESLineSearchMonitorSolutionUpdate", SNESLineSearchMonitorSolutionUpdate, NULL));
846bf7f4e0aSPeter Brune 
8471a9b3a06SPeter Brune   /* tolerances */
8489566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_minlambda", "Minimum step length", "SNESLineSearchSetTolerances", linesearch->steptol, &linesearch->steptol, NULL));
8499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_maxstep", "Maximum step size", "SNESLineSearchSetTolerances", linesearch->maxstep, &linesearch->maxstep, NULL));
8509566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_rtol", "Relative tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->rtol, &linesearch->rtol, NULL));
8519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_atol", "Absolute tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->atol, &linesearch->atol, NULL));
8529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_ltol", "Change in lambda tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->ltol, &linesearch->ltol, NULL));
8539566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_linesearch_max_it", "Maximum iterations for iterative line searches", "SNESLineSearchSetTolerances", linesearch->max_its, &linesearch->max_its, NULL));
8547a35526eSPeter Brune 
8551a9b3a06SPeter Brune   /* damping parameters */
8569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_damping", "Line search damping and initial step guess", "SNESLineSearchSetDamping", linesearch->damping, &linesearch->damping, NULL));
8571a9b3a06SPeter Brune 
8589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_linesearch_keeplambda", "Use previous lambda as damping", "SNESLineSearchSetKeepLambda", linesearch->keeplambda, &linesearch->keeplambda, NULL));
8591a9b3a06SPeter Brune 
8601a9b3a06SPeter Brune   /* precheck */
8619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_linesearch_precheck_picard", "Use a correction that sometimes improves convergence of Picard iteration", "SNESLineSearchPreCheckPicard", flg, &flg, &set));
8627a35526eSPeter Brune   if (set) {
8637a35526eSPeter Brune     if (flg) {
8647a35526eSPeter Brune       linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */
865f5af7f23SKarl Rupp 
866d0609cedSBarry 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));
8679566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchSetPreCheck(linesearch, SNESLineSearchPreCheckPicard, &linesearch->precheck_picard_angle));
8687a35526eSPeter Brune     } else {
8699566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchSetPreCheck(linesearch, NULL, NULL));
8707a35526eSPeter Brune     }
8717a35526eSPeter Brune   }
8729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_linesearch_order", "Order of approximation used in the line search", "SNESLineSearchSetOrder", linesearch->order, &linesearch->order, NULL));
8739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_linesearch_norms", "Compute final norms in line search", "SNESLineSearchSetComputeNorms", linesearch->norms, &linesearch->norms, NULL));
8747a35526eSPeter Brune 
875dbbe0bcdSBarry Smith   PetscTryTypeMethod(linesearch, setfromoptions, PetscOptionsObject);
8765a9b6599SPeter Brune 
877dbbe0bcdSBarry Smith   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)linesearch, PetscOptionsObject));
878d0609cedSBarry Smith   PetscOptionsEnd();
8793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
880bf7f4e0aSPeter Brune }
881bf7f4e0aSPeter Brune 
882f40b411bSPeter Brune /*@
883f190f2fcSBarry Smith   SNESLineSearchView - Prints useful information about the line search
884f40b411bSPeter Brune 
88520f4b53cSBarry Smith   Logically Collective
88620f4b53cSBarry Smith 
887f40b411bSPeter Brune   Input Parameters:
8882fe279fdSBarry Smith + linesearch - line search context
889*420bcc1bSBarry Smith - viewer     - the `PetscViewer` to display the line search information to
890f40b411bSPeter Brune 
891f40b411bSPeter Brune   Level: intermediate
892f40b411bSPeter Brune 
893*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `PetscViewer`, `SNESLineSearchCreate()`
894f40b411bSPeter Brune @*/
895d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer)
896d71ae5a4SJacob Faibussowitsch {
8977f1410a3SPeter Brune   PetscBool iascii;
898bf388a1fSBarry Smith 
899bf7f4e0aSPeter Brune   PetscFunctionBegin;
9007f1410a3SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
90148a46eb9SPierre Jolivet   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch), &viewer));
9027f1410a3SPeter Brune   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
9037f1410a3SPeter Brune   PetscCheckSameComm(linesearch, 1, viewer, 2);
904f40b411bSPeter Brune 
9059566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
9067f1410a3SPeter Brune   if (iascii) {
9079566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)linesearch, viewer));
9089566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
909dbbe0bcdSBarry Smith     PetscTryTypeMethod(linesearch, view, viewer);
9109566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
9119566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep, (double)linesearch->steptol));
9129566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol, (double)linesearch->atol, (double)linesearch->ltol));
91363a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  maximum iterations=%" PetscInt_FMT "\n", linesearch->max_its));
9146b2b7091SBarry Smith     if (linesearch->ops->precheck) {
9156b2b7091SBarry Smith       if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) {
91663a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  using precheck step to speed up Picard convergence\n"));
9177f1410a3SPeter Brune       } else {
91863a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  using user-defined precheck step\n"));
9197f1410a3SPeter Brune       }
9207f1410a3SPeter Brune     }
92148a46eb9SPierre Jolivet     if (linesearch->ops->postcheck) PetscCall(PetscViewerASCIIPrintf(viewer, "  using user-defined postcheck step\n"));
9227f1410a3SPeter Brune   }
9233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
924bf7f4e0aSPeter Brune }
925bf7f4e0aSPeter Brune 
926ea5d4fccSPeter Brune /*@C
927*420bcc1bSBarry Smith   SNESLineSearchGetType - Gets the `SNESLinesearchType` of a `SNESLineSearch`
928a80ff896SJed Brown 
929c3339decSBarry Smith   Logically Collective
930a80ff896SJed Brown 
9312fe279fdSBarry Smith   Input Parameter:
932*420bcc1bSBarry Smith . linesearch - the line search context
933a80ff896SJed Brown 
9342fe279fdSBarry Smith   Output Parameter:
9352fe279fdSBarry Smith . type - The type of line search, or `NULL` if not set
936a80ff896SJed Brown 
937a80ff896SJed Brown   Level: intermediate
938a80ff896SJed Brown 
939*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchSetType()`
940a80ff896SJed Brown @*/
941d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetType(SNESLineSearch linesearch, SNESLineSearchType *type)
942d71ae5a4SJacob Faibussowitsch {
943a80ff896SJed Brown   PetscFunctionBegin;
944a80ff896SJed Brown   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
9454f572ea9SToby Isaac   PetscAssertPointer(type, 2);
946a80ff896SJed Brown   *type = ((PetscObject)linesearch)->type_name;
9473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
948a80ff896SJed Brown }
949a80ff896SJed Brown 
950a80ff896SJed Brown /*@C
951*420bcc1bSBarry Smith   SNESLineSearchSetType - Sets the `SNESLinesearchType` of a `SNESLineSearch`
952f40b411bSPeter Brune 
953c3339decSBarry Smith   Logically Collective
954f190f2fcSBarry Smith 
955f40b411bSPeter Brune   Input Parameters:
956*420bcc1bSBarry Smith + linesearch - the line search context
957ceaaa498SBarry Smith - type       - The type of line search to be used, see `SNESLineSearchType`
9581fe24845SBarry Smith 
9593c7db156SBarry Smith   Options Database Key:
9600b00b554SBarry Smith . -snes_linesearch_type <type> - basic (or equivalently none), bt, l2, cp, nleqerr, shell
961cd7522eaSPeter Brune 
962f40b411bSPeter Brune   Level: intermediate
963f40b411bSPeter Brune 
964*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchGetType()`
965f40b411bSPeter Brune @*/
966d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type)
967d71ae5a4SJacob Faibussowitsch {
968bf7f4e0aSPeter Brune   PetscBool match;
9695f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(SNESLineSearch);
970bf7f4e0aSPeter Brune 
971bf7f4e0aSPeter Brune   PetscFunctionBegin;
972f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
9734f572ea9SToby Isaac   PetscAssertPointer(type, 2);
974bf7f4e0aSPeter Brune 
9759566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)linesearch, type, &match));
9763ba16761SJacob Faibussowitsch   if (match) PetscFunctionReturn(PETSC_SUCCESS);
977bf7f4e0aSPeter Brune 
9789566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(SNESLineSearchList, type, &r));
9796adde796SStefano Zampini   PetscCheck(r, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested Line Search type %s", type);
980bf7f4e0aSPeter Brune   /* Destroy the previous private line search context */
981bf7f4e0aSPeter Brune   if (linesearch->ops->destroy) {
9829566063dSJacob Faibussowitsch     PetscCall((*(linesearch)->ops->destroy)(linesearch));
9830298fd71SBarry Smith     linesearch->ops->destroy = NULL;
984bf7f4e0aSPeter Brune   }
985f1c6b773SPeter Brune   /* Reinitialize function pointers in SNESLineSearchOps structure */
9869e5d0892SLisandro Dalcin   linesearch->ops->apply          = NULL;
9879e5d0892SLisandro Dalcin   linesearch->ops->view           = NULL;
9889e5d0892SLisandro Dalcin   linesearch->ops->setfromoptions = NULL;
9899e5d0892SLisandro Dalcin   linesearch->ops->destroy        = NULL;
990bf7f4e0aSPeter Brune 
9919566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)linesearch, type));
9929566063dSJacob Faibussowitsch   PetscCall((*r)(linesearch));
9933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
994bf7f4e0aSPeter Brune }
995bf7f4e0aSPeter Brune 
996f40b411bSPeter Brune /*@
997f6dfbefdSBarry Smith   SNESLineSearchSetSNES - Sets the `SNES` for the linesearch for function evaluation.
998f40b411bSPeter Brune 
999f40b411bSPeter Brune   Input Parameters:
1000*420bcc1bSBarry Smith + linesearch - the line search context
1001*420bcc1bSBarry Smith - snes       - The `SNES` instance
1002f40b411bSPeter Brune 
100378bcb3b5SPeter Brune   Level: developer
100478bcb3b5SPeter Brune 
1005f6dfbefdSBarry Smith   Note:
1006f190f2fcSBarry Smith   This happens automatically when the line search is obtained/created with
1007f6dfbefdSBarry Smith   `SNESGetLineSearch()`.  This routine is therefore mainly called within `SNES`
100878bcb3b5SPeter Brune   implementations.
1009f40b411bSPeter Brune 
1010*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetSNES()`, `SNESLineSearchSetVecs()`
1011f40b411bSPeter Brune @*/
1012d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes)
1013d71ae5a4SJacob Faibussowitsch {
1014bf7f4e0aSPeter Brune   PetscFunctionBegin;
1015f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1016bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 2);
1017bf7f4e0aSPeter Brune   linesearch->snes = snes;
10183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1019bf7f4e0aSPeter Brune }
1020bf7f4e0aSPeter Brune 
1021f40b411bSPeter Brune /*@
1022f6dfbefdSBarry Smith   SNESLineSearchGetSNES - Gets the `SNES` instance associated with the line search.
1023f6dfbefdSBarry Smith 
1024f6dfbefdSBarry Smith   Not Collective
1025f40b411bSPeter Brune 
10262fe279fdSBarry Smith   Input Parameter:
1027*420bcc1bSBarry Smith . linesearch - the line search context
1028f40b411bSPeter Brune 
10292fe279fdSBarry Smith   Output Parameter:
10302fe279fdSBarry Smith . snes - The `SNES` instance
1031f40b411bSPeter Brune 
10328141a3b9SPeter Brune   Level: developer
1033f40b411bSPeter Brune 
1034*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESType`, `SNESLineSearchSetVecs()`
1035f40b411bSPeter Brune @*/
1036d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes)
1037d71ae5a4SJacob Faibussowitsch {
1038bf7f4e0aSPeter Brune   PetscFunctionBegin;
1039f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
10404f572ea9SToby Isaac   PetscAssertPointer(snes, 2);
1041bf7f4e0aSPeter Brune   *snes = linesearch->snes;
10423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1043bf7f4e0aSPeter Brune }
1044bf7f4e0aSPeter Brune 
1045f40b411bSPeter Brune /*@
1046*420bcc1bSBarry Smith   SNESLineSearchGetLambda - Gets the last line search steplength used
1047f40b411bSPeter Brune 
1048f6dfbefdSBarry Smith   Not Collective
1049f6dfbefdSBarry Smith 
10502fe279fdSBarry Smith   Input Parameter:
1051*420bcc1bSBarry Smith . linesearch - the line search context
1052f40b411bSPeter Brune 
10532fe279fdSBarry Smith   Output Parameter:
1054f6dfbefdSBarry Smith . lambda - The last steplength computed during `SNESLineSearchApply()`
1055f40b411bSPeter Brune 
105678bcb3b5SPeter Brune   Level: advanced
105778bcb3b5SPeter Brune 
1058f6dfbefdSBarry Smith   Note:
10598e557f58SPeter Brune   This is useful in methods where the solver is ill-scaled and
106078bcb3b5SPeter Brune   requires some adaptive notion of the difference in scale between the
1061f6dfbefdSBarry Smith   solution and the function.  For instance, `SNESQN` may be scaled by the
106278bcb3b5SPeter Brune   line search lambda using the argument -snes_qn_scaling ls.
106378bcb3b5SPeter Brune 
1064*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetLambda()`, `SNESLineSearchGetDamping()`, `SNESLineSearchApply()`
1065f40b411bSPeter Brune @*/
1066d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch linesearch, PetscReal *lambda)
1067d71ae5a4SJacob Faibussowitsch {
10686a388c36SPeter Brune   PetscFunctionBegin;
1069f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
10704f572ea9SToby Isaac   PetscAssertPointer(lambda, 2);
10716a388c36SPeter Brune   *lambda = linesearch->lambda;
10723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10736a388c36SPeter Brune }
10746a388c36SPeter Brune 
1075f40b411bSPeter Brune /*@
1076f6dfbefdSBarry Smith   SNESLineSearchSetLambda - Sets the line search steplength
1077f40b411bSPeter Brune 
1078f40b411bSPeter Brune   Input Parameters:
10798e557f58SPeter Brune + linesearch - line search context
1080*420bcc1bSBarry Smith - lambda     - The steplength to use
1081f40b411bSPeter Brune 
108220f4b53cSBarry Smith   Level: advanced
108320f4b53cSBarry Smith 
1084f6dfbefdSBarry Smith   Note:
1085f6dfbefdSBarry Smith   This routine is typically used within implementations of `SNESLineSearchApply()`
1086f6dfbefdSBarry Smith   to set the final steplength.  This routine (and `SNESLineSearchGetLambda()`) were
1087cd7522eaSPeter Brune   added in order to facilitate Quasi-Newton methods that use the previous steplength
1088cd7522eaSPeter Brune   as an inner scaling parameter.
1089cd7522eaSPeter Brune 
1090*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetLambda()`
1091f40b411bSPeter Brune @*/
1092d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
1093d71ae5a4SJacob Faibussowitsch {
10946a388c36SPeter Brune   PetscFunctionBegin;
1095f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
10966a388c36SPeter Brune   linesearch->lambda = lambda;
10973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10986a388c36SPeter Brune }
10996a388c36SPeter Brune 
1100f40b411bSPeter Brune /*@
1101ceaaa498SBarry Smith   SNESLineSearchGetTolerances - Gets the tolerances for the line search.
1102f40b411bSPeter Brune 
1103f6dfbefdSBarry Smith   Not Collective
1104f6dfbefdSBarry Smith 
1105f899ff85SJose E. Roman   Input Parameter:
1106*420bcc1bSBarry Smith . linesearch - the line search context
1107f40b411bSPeter Brune 
1108f40b411bSPeter Brune   Output Parameters:
1109b13b64c2SBarry Smith + steptol - The minimum steplength
11106cc8e53bSPeter Brune . maxstep - The maximum steplength
1111516fe3c3SPeter Brune . rtol    - The relative tolerance for iterative line searches
1112516fe3c3SPeter Brune . atol    - The absolute tolerance for iterative line searches
1113516fe3c3SPeter Brune . ltol    - The change in lambda tolerance for iterative line searches
1114e4094ef1SJacob Faibussowitsch - max_its - The maximum number of iterations of the line search
1115f40b411bSPeter Brune 
111678bcb3b5SPeter Brune   Level: intermediate
111778bcb3b5SPeter Brune 
1118f6dfbefdSBarry Smith   Note:
111978bcb3b5SPeter Brune   Different line searches may implement these parameters slightly differently as
11203c7d6663SPeter Brune   the type requires.
1121516fe3c3SPeter Brune 
1122*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetTolerances()`
1123f40b411bSPeter Brune @*/
1124d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch linesearch, PetscReal *steptol, PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
1125d71ae5a4SJacob Faibussowitsch {
11266a388c36SPeter Brune   PetscFunctionBegin;
1127f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1128516fe3c3SPeter Brune   if (steptol) {
11294f572ea9SToby Isaac     PetscAssertPointer(steptol, 2);
11306a388c36SPeter Brune     *steptol = linesearch->steptol;
1131516fe3c3SPeter Brune   }
1132516fe3c3SPeter Brune   if (maxstep) {
11334f572ea9SToby Isaac     PetscAssertPointer(maxstep, 3);
1134b13b64c2SBarry Smith     *maxstep = linesearch->maxstep;
1135516fe3c3SPeter Brune   }
1136516fe3c3SPeter Brune   if (rtol) {
11374f572ea9SToby Isaac     PetscAssertPointer(rtol, 4);
1138516fe3c3SPeter Brune     *rtol = linesearch->rtol;
1139516fe3c3SPeter Brune   }
1140516fe3c3SPeter Brune   if (atol) {
11414f572ea9SToby Isaac     PetscAssertPointer(atol, 5);
1142516fe3c3SPeter Brune     *atol = linesearch->atol;
1143516fe3c3SPeter Brune   }
1144516fe3c3SPeter Brune   if (ltol) {
11454f572ea9SToby Isaac     PetscAssertPointer(ltol, 6);
1146516fe3c3SPeter Brune     *ltol = linesearch->ltol;
1147516fe3c3SPeter Brune   }
1148516fe3c3SPeter Brune   if (max_its) {
11494f572ea9SToby Isaac     PetscAssertPointer(max_its, 7);
1150516fe3c3SPeter Brune     *max_its = linesearch->max_its;
1151516fe3c3SPeter Brune   }
11523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11536a388c36SPeter Brune }
11546a388c36SPeter Brune 
1155f40b411bSPeter Brune /*@
1156ceaaa498SBarry Smith   SNESLineSearchSetTolerances -  Sets the tolerances for the linesearch.
1157f40b411bSPeter Brune 
1158f6dfbefdSBarry Smith   Collective
1159f6dfbefdSBarry Smith 
1160f40b411bSPeter Brune   Input Parameters:
1161*420bcc1bSBarry Smith + linesearch - the line search context
1162516fe3c3SPeter Brune . steptol    - The minimum steplength
11636cc8e53bSPeter Brune . maxstep    - The maximum steplength
1164516fe3c3SPeter Brune . rtol       - The relative tolerance for iterative line searches
1165516fe3c3SPeter Brune . atol       - The absolute tolerance for iterative line searches
1166516fe3c3SPeter Brune . ltol       - The change in lambda tolerance for iterative line searches
1167*420bcc1bSBarry Smith - max_it     - The maximum number of iterations of the line search
1168*420bcc1bSBarry Smith 
1169*420bcc1bSBarry Smith   Options Database Keys:
1170*420bcc1bSBarry Smith + -snes_linesearch_minlambda - The minimum step length
1171*420bcc1bSBarry Smith . -snes_linesearch_maxstep   - The maximum step size
1172*420bcc1bSBarry Smith . -snes_linesearch_rtol      - Relative tolerance for iterative line searches
1173*420bcc1bSBarry Smith . -snes_linesearch_atol      - Absolute tolerance for iterative line searches
1174*420bcc1bSBarry Smith . -snes_linesearch_ltol      - Change in lambda tolerance for iterative line searches
1175*420bcc1bSBarry Smith - -snes_linesearch_max_it    - The number of iterations for iterative line searches
1176f40b411bSPeter Brune 
117720f4b53cSBarry Smith   Level: intermediate
117820f4b53cSBarry Smith 
1179f6dfbefdSBarry Smith   Note:
1180*420bcc1bSBarry Smith   The user may choose to not set any of the tolerances using `PETSC_DEFAULT` in place of an argument.
1181f40b411bSPeter Brune 
1182*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetTolerances()`
1183f40b411bSPeter Brune @*/
1184*420bcc1bSBarry Smith PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch linesearch, PetscReal steptol, PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_it)
1185d71ae5a4SJacob Faibussowitsch {
11866a388c36SPeter Brune   PetscFunctionBegin;
1187f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1188d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, steptol, 2);
1189d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, maxstep, 3);
1190d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, rtol, 4);
1191d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, atol, 5);
1192d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, ltol, 6);
1193*420bcc1bSBarry Smith   PetscValidLogicalCollectiveInt(linesearch, max_it, 7);
1194d3952184SSatish Balay 
119513bcc0bdSJacob Faibussowitsch   if (steptol != (PetscReal)PETSC_DEFAULT) {
11965f80ce2aSJacob Faibussowitsch     PetscCheck(steptol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Minimum step length %14.12e must be non-negative", (double)steptol);
11976a388c36SPeter Brune     linesearch->steptol = steptol;
1198d3952184SSatish Balay   }
1199d3952184SSatish Balay 
120013bcc0bdSJacob Faibussowitsch   if (maxstep != (PetscReal)PETSC_DEFAULT) {
12015f80ce2aSJacob Faibussowitsch     PetscCheck(maxstep >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum step length %14.12e must be non-negative", (double)maxstep);
1202516fe3c3SPeter Brune     linesearch->maxstep = maxstep;
1203d3952184SSatish Balay   }
1204d3952184SSatish Balay 
120513bcc0bdSJacob Faibussowitsch   if (rtol != (PetscReal)PETSC_DEFAULT) {
12062061ca32SJunchao 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);
1207516fe3c3SPeter Brune     linesearch->rtol = rtol;
1208d3952184SSatish Balay   }
1209d3952184SSatish Balay 
121013bcc0bdSJacob Faibussowitsch   if (atol != (PetscReal)PETSC_DEFAULT) {
12115f80ce2aSJacob Faibussowitsch     PetscCheck(atol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Absolute tolerance %14.12e must be non-negative", (double)atol);
1212516fe3c3SPeter Brune     linesearch->atol = atol;
1213d3952184SSatish Balay   }
1214d3952184SSatish Balay 
121513bcc0bdSJacob Faibussowitsch   if (ltol != (PetscReal)PETSC_DEFAULT) {
12165f80ce2aSJacob Faibussowitsch     PetscCheck(ltol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Lambda tolerance %14.12e must be non-negative", (double)ltol);
1217516fe3c3SPeter Brune     linesearch->ltol = ltol;
1218d3952184SSatish Balay   }
1219d3952184SSatish Balay 
1220*420bcc1bSBarry Smith   if (max_it != PETSC_DEFAULT) {
1221*420bcc1bSBarry Smith     PetscCheck(max_it >= 0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of iterations %" PetscInt_FMT " must be non-negative", max_it);
1222*420bcc1bSBarry Smith     linesearch->max_its = max_it;
1223d3952184SSatish Balay   }
12243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12256a388c36SPeter Brune }
12266a388c36SPeter Brune 
1227f40b411bSPeter Brune /*@
1228f1c6b773SPeter Brune   SNESLineSearchGetDamping - Gets the line search damping parameter.
1229f40b411bSPeter Brune 
12302fe279fdSBarry Smith   Input Parameter:
1231*420bcc1bSBarry Smith . linesearch - the line search context
1232f40b411bSPeter Brune 
12332fe279fdSBarry Smith   Output Parameter:
12348e557f58SPeter Brune . damping - The damping parameter
1235f40b411bSPeter Brune 
123678bcb3b5SPeter Brune   Level: advanced
1237f40b411bSPeter Brune 
1238*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearchGetStepTolerance()`, `SNESQN`
1239f40b411bSPeter Brune @*/
1240d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch linesearch, PetscReal *damping)
1241d71ae5a4SJacob Faibussowitsch {
12426a388c36SPeter Brune   PetscFunctionBegin;
1243f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
12444f572ea9SToby Isaac   PetscAssertPointer(damping, 2);
12456a388c36SPeter Brune   *damping = linesearch->damping;
12463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12476a388c36SPeter Brune }
12486a388c36SPeter Brune 
1249f40b411bSPeter Brune /*@
1250fd292e60Sprj-   SNESLineSearchSetDamping - Sets the line search damping parameter.
1251f40b411bSPeter Brune 
1252f40b411bSPeter Brune   Input Parameters:
1253*420bcc1bSBarry Smith + linesearch - the line search context
125403fd4120SBarry Smith - damping    - The damping parameter
1255f40b411bSPeter Brune 
1256f6dfbefdSBarry Smith   Options Database Key:
1257f6dfbefdSBarry Smith . -snes_linesearch_damping <damping> - the damping value
1258f6dfbefdSBarry Smith 
1259f40b411bSPeter Brune   Level: intermediate
1260f40b411bSPeter Brune 
1261f6dfbefdSBarry Smith   Note:
1262f6dfbefdSBarry Smith   The `SNESLINESEARCHNONE` line search merely takes the update step scaled by the damping parameter.
1263*420bcc1bSBarry Smith   The use of the damping parameter in the `SNESLINESEARCHL2` and `SNESLINESEARCHCP` line searches is much more subtle;
126478bcb3b5SPeter Brune   it is used as a starting point in calculating the secant step. However, the eventual
1265*420bcc1bSBarry Smith   step may be of greater length than the damping parameter.  In the `SNESLINESEARCHBT` line search it is
1266*420bcc1bSBarry Smith   used as the maximum possible step length, as the `SNESLINESEARCHBT` line search only backtracks.
1267cd7522eaSPeter Brune 
1268*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetDamping()`
1269f40b411bSPeter Brune @*/
1270d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch linesearch, PetscReal damping)
1271d71ae5a4SJacob Faibussowitsch {
12726a388c36SPeter Brune   PetscFunctionBegin;
1273f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
12746a388c36SPeter Brune   linesearch->damping = damping;
12753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12766a388c36SPeter Brune }
12776a388c36SPeter Brune 
127859405d5eSPeter Brune /*@
127959405d5eSPeter Brune   SNESLineSearchGetOrder - Gets the line search approximation order.
128059405d5eSPeter Brune 
1281f6dfbefdSBarry Smith   Input Parameter:
1282*420bcc1bSBarry Smith . linesearch - the line search context
128359405d5eSPeter Brune 
1284f6dfbefdSBarry Smith   Output Parameter:
12858e557f58SPeter Brune . order - The order
128659405d5eSPeter Brune 
128759405d5eSPeter Brune   Level: intermediate
128859405d5eSPeter Brune 
1289*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetOrder()`
129059405d5eSPeter Brune @*/
1291d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch linesearch, PetscInt *order)
1292d71ae5a4SJacob Faibussowitsch {
129359405d5eSPeter Brune   PetscFunctionBegin;
129459405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
12954f572ea9SToby Isaac   PetscAssertPointer(order, 2);
129659405d5eSPeter Brune   *order = linesearch->order;
12973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
129859405d5eSPeter Brune }
129959405d5eSPeter Brune 
130059405d5eSPeter Brune /*@
13011f8196a2SJed Brown   SNESLineSearchSetOrder - Sets the maximum order of the polynomial fit used in the line search
130259405d5eSPeter Brune 
130359405d5eSPeter Brune   Input Parameters:
1304*420bcc1bSBarry Smith + linesearch - the line search context
1305ceaaa498SBarry Smith - order      - The order
130659405d5eSPeter Brune 
130759405d5eSPeter Brune   Level: intermediate
130859405d5eSPeter Brune 
1309ceaaa498SBarry Smith   Values for `order`\:
1310f6dfbefdSBarry Smith +  1 or `SNES_LINESEARCH_ORDER_LINEAR` - linear order
1311f6dfbefdSBarry Smith .  2 or `SNES_LINESEARCH_ORDER_QUADRATIC` - quadratic order
1312f6dfbefdSBarry Smith -  3 or `SNES_LINESEARCH_ORDER_CUBIC` - cubic order
131378bcb3b5SPeter Brune 
1314*420bcc1bSBarry Smith   Options Database Key:
1315*420bcc1bSBarry Smith . -snes_linesearch_order <order> - 1, 2, 3.  Most types only support certain orders (`SNESLINESEARCHBT` supports 2 or 3)
1316*420bcc1bSBarry Smith 
1317ceaaa498SBarry Smith   Note:
1318ceaaa498SBarry Smith   These orders are supported by `SNESLINESEARCHBT` and `SNESLINESEARCHCP`
131959405d5eSPeter Brune 
1320*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetOrder()`, `SNESLineSearchSetDamping()`
132159405d5eSPeter Brune @*/
1322d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch linesearch, PetscInt order)
1323d71ae5a4SJacob Faibussowitsch {
132459405d5eSPeter Brune   PetscFunctionBegin;
132559405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
132659405d5eSPeter Brune   linesearch->order = order;
13273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
132859405d5eSPeter Brune }
132959405d5eSPeter Brune 
1330f40b411bSPeter Brune /*@
1331*420bcc1bSBarry Smith   SNESLineSearchGetNorms - Gets the norms for the current solution `X`, the current update `Y`, and the current function value `F`.
1332f40b411bSPeter Brune 
1333f6dfbefdSBarry Smith   Not Collective
1334f6dfbefdSBarry Smith 
1335f899ff85SJose E. Roman   Input Parameter:
1336*420bcc1bSBarry Smith . linesearch - the line search context
1337f40b411bSPeter Brune 
1338f40b411bSPeter Brune   Output Parameters:
1339f40b411bSPeter Brune + xnorm - The norm of the current solution
1340a68bbae5SBarry Smith . fnorm - The norm of the current function, this is the `norm(function(X))` where `X` is the current solution.
1341f40b411bSPeter Brune - ynorm - The norm of the current update
1342f40b411bSPeter Brune 
134378bcb3b5SPeter Brune   Level: developer
1344f40b411bSPeter Brune 
1345*420bcc1bSBarry Smith   Notes:
1346*420bcc1bSBarry Smith   Some values may not be up-to-date at particular points in the code.
1347a68bbae5SBarry Smith 
1348a68bbae5SBarry Smith   This, in combination with `SNESLineSearchSetNorms()`, allow the line search and the `SNESSolve_XXX()` to share
1349a68bbae5SBarry Smith   computed values.
1350a68bbae5SBarry Smith 
1351*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetNorms()` `SNESLineSearchGetVecs()`
1352f40b411bSPeter Brune @*/
1353d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal *xnorm, PetscReal *fnorm, PetscReal *ynorm)
1354d71ae5a4SJacob Faibussowitsch {
1355bf7f4e0aSPeter Brune   PetscFunctionBegin;
1356f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1357f5af7f23SKarl Rupp   if (xnorm) *xnorm = linesearch->xnorm;
1358f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
1359f5af7f23SKarl Rupp   if (ynorm) *ynorm = linesearch->ynorm;
13603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1361bf7f4e0aSPeter Brune }
1362bf7f4e0aSPeter Brune 
1363f40b411bSPeter Brune /*@
1364*420bcc1bSBarry Smith   SNESLineSearchSetNorms - Gets the computed norms for the current solution `X`, the current update `Y`, and the current function value `F`.
1365f40b411bSPeter Brune 
1366c3339decSBarry Smith   Collective
1367f6dfbefdSBarry Smith 
1368f40b411bSPeter Brune   Input Parameters:
1369*420bcc1bSBarry Smith + linesearch - the line search context
1370f40b411bSPeter Brune . xnorm      - The norm of the current solution
1371a68bbae5SBarry Smith . fnorm      - The norm of the current function, this is the `norm(function(X))` where `X` is the current solution
1372f40b411bSPeter Brune - ynorm      - The norm of the current update
1373f40b411bSPeter Brune 
1374f6dfbefdSBarry Smith   Level: developer
1375f40b411bSPeter Brune 
1376*420bcc1bSBarry Smith   Note:
1377*420bcc1bSBarry Smith   This is called by the line search routines to store the values they have just computed
1378*420bcc1bSBarry Smith 
1379*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()`
1380f40b411bSPeter Brune @*/
1381d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
1382d71ae5a4SJacob Faibussowitsch {
13836a388c36SPeter Brune   PetscFunctionBegin;
1384f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
13856a388c36SPeter Brune   linesearch->xnorm = xnorm;
13866a388c36SPeter Brune   linesearch->fnorm = fnorm;
13876a388c36SPeter Brune   linesearch->ynorm = ynorm;
13883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13896a388c36SPeter Brune }
13906a388c36SPeter Brune 
1391f40b411bSPeter Brune /*@
1392*420bcc1bSBarry Smith   SNESLineSearchComputeNorms - Explicitly computes the norms of the current solution `X`, the current update `Y`, and the current function value `F`.
1393f40b411bSPeter Brune 
13942fe279fdSBarry Smith   Input Parameter:
1395*420bcc1bSBarry Smith . linesearch - the line search context
1396f40b411bSPeter Brune 
139720f4b53cSBarry Smith   Options Database Key:
13988e557f58SPeter Brune . -snes_linesearch_norms - turn norm computation on or off
1399f40b411bSPeter Brune 
1400f40b411bSPeter Brune   Level: intermediate
1401f40b411bSPeter Brune 
1402*420bcc1bSBarry Smith   Developer Note:
1403*420bcc1bSBarry Smith   The options database key is misnamed. It should be -snes_linesearch_compute_norms
1404*420bcc1bSBarry Smith 
1405*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms`, `SNESLineSearchSetNorms()`, `SNESLineSearchSetComputeNorms()`
1406f40b411bSPeter Brune @*/
1407d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
1408d71ae5a4SJacob Faibussowitsch {
14099bd66eb0SPeter Brune   SNES snes;
1410bf388a1fSBarry Smith 
14116a388c36SPeter Brune   PetscFunctionBegin;
14126a388c36SPeter Brune   if (linesearch->norms) {
14139bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
14149566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
14159566063dSJacob Faibussowitsch       PetscCall(VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
14169566063dSJacob Faibussowitsch       PetscCall(VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm));
14179566063dSJacob Faibussowitsch       PetscCall((*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm));
14189bd66eb0SPeter Brune     } else {
14199566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(linesearch->vec_func, NORM_2, &linesearch->fnorm));
14209566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
14219566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm));
14229566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(linesearch->vec_func, NORM_2, &linesearch->fnorm));
14239566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
14249566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(linesearch->vec_update, NORM_2, &linesearch->ynorm));
14256a388c36SPeter Brune     }
14269bd66eb0SPeter Brune   }
14273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14286a388c36SPeter Brune }
14296a388c36SPeter Brune 
14306f263ca3SPeter Brune /*@
14316f263ca3SPeter Brune   SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search.
14326f263ca3SPeter Brune 
14336f263ca3SPeter Brune   Input Parameters:
1434*420bcc1bSBarry Smith + linesearch - the line search context
143578bcb3b5SPeter Brune - flg        - indicates whether or not to compute norms
14366f263ca3SPeter Brune 
143720f4b53cSBarry Smith   Options Database Key:
1438*420bcc1bSBarry Smith . -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic (none) `SNESLINESEARCHBASIC` line search
14396f263ca3SPeter Brune 
144020f4b53cSBarry Smith   Level: intermediate
144120f4b53cSBarry Smith 
1442f6dfbefdSBarry Smith   Note:
1443f6dfbefdSBarry 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.
14446f263ca3SPeter Brune 
1445*420bcc1bSBarry Smith   Developer Note:
1446*420bcc1bSBarry Smith   The options database key is misnamed. It should be -snes_linesearch_compute_norms
1447*420bcc1bSBarry Smith 
1448*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetNorms()`, `SNESLineSearchComputeNorms()`, `SNESLINESEARCHBASIC`
14496f263ca3SPeter Brune @*/
1450d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg)
1451d71ae5a4SJacob Faibussowitsch {
14526f263ca3SPeter Brune   PetscFunctionBegin;
14536f263ca3SPeter Brune   linesearch->norms = flg;
14543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14556f263ca3SPeter Brune }
14566f263ca3SPeter Brune 
1457f40b411bSPeter Brune /*@
1458f6dfbefdSBarry Smith   SNESLineSearchGetVecs - Gets the vectors from the `SNESLineSearch` context
1459f6dfbefdSBarry Smith 
14608f14a041SBarry Smith   Not Collective but the vectors are parallel
1461f40b411bSPeter Brune 
1462f899ff85SJose E. Roman   Input Parameter:
1463*420bcc1bSBarry Smith . linesearch - the line search context
1464f40b411bSPeter Brune 
1465f40b411bSPeter Brune   Output Parameters:
14666232e825SPeter Brune + X - Solution vector
14676232e825SPeter Brune . F - Function vector
14686232e825SPeter Brune . Y - Search direction vector
14696232e825SPeter Brune . W - Solution work vector
14706232e825SPeter Brune - G - Function work vector
14716232e825SPeter Brune 
147220f4b53cSBarry Smith   Level: advanced
147320f4b53cSBarry Smith 
14747bba9028SPeter Brune   Notes:
147520f4b53cSBarry Smith   At the beginning of a line search application, `X` should contain a
147620f4b53cSBarry Smith   solution and the vector `F` the function computed at `X`.  At the end of the
147720f4b53cSBarry Smith   line search application, `X` should contain the new solution, and `F` the
14786232e825SPeter Brune   function evaluated at the new solution.
1479f40b411bSPeter Brune 
1480f6dfbefdSBarry Smith   These vectors are owned by the `SNESLineSearch` and should not be destroyed by the caller
14812a7a6963SBarry Smith 
1482*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()`
1483f40b411bSPeter Brune @*/
1484d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch, Vec *X, Vec *F, Vec *Y, Vec *W, Vec *G)
1485d71ae5a4SJacob Faibussowitsch {
14866a388c36SPeter Brune   PetscFunctionBegin;
1487f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
14886a388c36SPeter Brune   if (X) {
14894f572ea9SToby Isaac     PetscAssertPointer(X, 2);
14906a388c36SPeter Brune     *X = linesearch->vec_sol;
14916a388c36SPeter Brune   }
14926a388c36SPeter Brune   if (F) {
14934f572ea9SToby Isaac     PetscAssertPointer(F, 3);
14946a388c36SPeter Brune     *F = linesearch->vec_func;
14956a388c36SPeter Brune   }
14966a388c36SPeter Brune   if (Y) {
14974f572ea9SToby Isaac     PetscAssertPointer(Y, 4);
14986a388c36SPeter Brune     *Y = linesearch->vec_update;
14996a388c36SPeter Brune   }
15006a388c36SPeter Brune   if (W) {
15014f572ea9SToby Isaac     PetscAssertPointer(W, 5);
15026a388c36SPeter Brune     *W = linesearch->vec_sol_new;
15036a388c36SPeter Brune   }
15046a388c36SPeter Brune   if (G) {
15054f572ea9SToby Isaac     PetscAssertPointer(G, 6);
15066a388c36SPeter Brune     *G = linesearch->vec_func_new;
15076a388c36SPeter Brune   }
15083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15096a388c36SPeter Brune }
15106a388c36SPeter Brune 
1511f40b411bSPeter Brune /*@
1512f6dfbefdSBarry Smith   SNESLineSearchSetVecs - Sets the vectors on the `SNESLineSearch` context
1513f6dfbefdSBarry Smith 
1514c3339decSBarry Smith   Logically Collective
1515f40b411bSPeter Brune 
1516f40b411bSPeter Brune   Input Parameters:
1517*420bcc1bSBarry Smith + linesearch - the line search context
15186232e825SPeter Brune . X          - Solution vector
15196232e825SPeter Brune . F          - Function vector
15206232e825SPeter Brune . Y          - Search direction vector
15216232e825SPeter Brune . W          - Solution work vector
15226232e825SPeter Brune - G          - Function work vector
1523f40b411bSPeter Brune 
1524*420bcc1bSBarry Smith   Level: developer
1525f40b411bSPeter Brune 
1526*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetNorms()`, `SNESLineSearchGetVecs()`
1527f40b411bSPeter Brune @*/
1528d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch, Vec X, Vec F, Vec Y, Vec W, Vec G)
1529d71ae5a4SJacob Faibussowitsch {
15306a388c36SPeter Brune   PetscFunctionBegin;
1531f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
15326a388c36SPeter Brune   if (X) {
15336a388c36SPeter Brune     PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
15346a388c36SPeter Brune     linesearch->vec_sol = X;
15356a388c36SPeter Brune   }
15366a388c36SPeter Brune   if (F) {
15376a388c36SPeter Brune     PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
15386a388c36SPeter Brune     linesearch->vec_func = F;
15396a388c36SPeter Brune   }
15406a388c36SPeter Brune   if (Y) {
15416a388c36SPeter Brune     PetscValidHeaderSpecific(Y, VEC_CLASSID, 4);
15426a388c36SPeter Brune     linesearch->vec_update = Y;
15436a388c36SPeter Brune   }
15446a388c36SPeter Brune   if (W) {
15456a388c36SPeter Brune     PetscValidHeaderSpecific(W, VEC_CLASSID, 5);
15466a388c36SPeter Brune     linesearch->vec_sol_new = W;
15476a388c36SPeter Brune   }
15486a388c36SPeter Brune   if (G) {
15496a388c36SPeter Brune     PetscValidHeaderSpecific(G, VEC_CLASSID, 6);
15506a388c36SPeter Brune     linesearch->vec_func_new = G;
15516a388c36SPeter Brune   }
15523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15536a388c36SPeter Brune }
15546a388c36SPeter Brune 
1555e7058c64SPeter Brune /*@C
1556f1c6b773SPeter Brune   SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1557f6dfbefdSBarry Smith   `SNESLineSearch` options in the database.
1558e7058c64SPeter Brune 
1559c3339decSBarry Smith   Logically Collective
1560e7058c64SPeter Brune 
1561e7058c64SPeter Brune   Input Parameters:
1562f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context
1563e7058c64SPeter Brune - prefix     - the prefix to prepend to all option names
1564e7058c64SPeter Brune 
156520f4b53cSBarry Smith   Level: advanced
156620f4b53cSBarry Smith 
1567f6dfbefdSBarry Smith   Note:
1568e7058c64SPeter Brune   A hyphen (-) must NOT be given at the beginning of the prefix name.
1569e7058c64SPeter Brune   The first character of all runtime options is AUTOMATICALLY the hyphen.
1570e7058c64SPeter Brune 
1571*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch()`, `SNESLineSearchSetFromOptions()`, `SNESGetOptionsPrefix()`
1572e7058c64SPeter Brune @*/
1573d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch, const char prefix[])
1574d71ae5a4SJacob Faibussowitsch {
1575e7058c64SPeter Brune   PetscFunctionBegin;
1576f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
15779566063dSJacob Faibussowitsch   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)linesearch, prefix));
15783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1579e7058c64SPeter Brune }
1580e7058c64SPeter Brune 
1581e7058c64SPeter Brune /*@C
1582f6dfbefdSBarry Smith   SNESLineSearchGetOptionsPrefix - Gets the prefix used for searching for all
1583f1c6b773SPeter Brune   SNESLineSearch options in the database.
1584e7058c64SPeter Brune 
1585e7058c64SPeter Brune   Not Collective
1586e7058c64SPeter Brune 
1587e7058c64SPeter Brune   Input Parameter:
1588f6dfbefdSBarry Smith . linesearch - the `SNESLineSearch` context
1589e7058c64SPeter Brune 
1590e7058c64SPeter Brune   Output Parameter:
1591e7058c64SPeter Brune . prefix - pointer to the prefix string used
1592e7058c64SPeter Brune 
1593e7058c64SPeter Brune   Level: advanced
1594e7058c64SPeter Brune 
1595e4094ef1SJacob Faibussowitsch   Fortran Notes:
159620f4b53cSBarry Smith   The user should pass in a string 'prefix' of
159720f4b53cSBarry Smith   sufficient length to hold the prefix.
159820f4b53cSBarry Smith 
1599*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESAppendOptionsPrefix()`
1600e7058c64SPeter Brune @*/
1601d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch, const char *prefix[])
1602d71ae5a4SJacob Faibussowitsch {
1603e7058c64SPeter Brune   PetscFunctionBegin;
1604f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
16059566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)linesearch, prefix));
16063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1607e7058c64SPeter Brune }
1608bf7f4e0aSPeter Brune 
16098d359177SBarry Smith /*@C
1610f6dfbefdSBarry Smith   SNESLineSearchSetWorkVecs - Sets work vectors for the line search.
1611f40b411bSPeter Brune 
1612d8d19677SJose E. Roman   Input Parameters:
1613f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context
1614f40b411bSPeter Brune - nwork      - the number of work vectors
1615f40b411bSPeter Brune 
1616f40b411bSPeter Brune   Level: developer
1617f40b411bSPeter Brune 
1618*420bcc1bSBarry Smith   Developer Note:
1619*420bcc1bSBarry Smith   This is called from within the set up routines for each of the line search types `SNESLineSearchType`
1620*420bcc1bSBarry Smith 
1621*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESSetWorkVecs()`
1622f40b411bSPeter Brune @*/
1623d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork)
1624d71ae5a4SJacob Faibussowitsch {
1625bf7f4e0aSPeter Brune   PetscFunctionBegin;
16260fdf79fbSJacob Faibussowitsch   PetscCheck(linesearch->vec_sol, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
16279566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work));
16283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1629bf7f4e0aSPeter Brune }
1630bf7f4e0aSPeter Brune 
1631f40b411bSPeter Brune /*@
1632422a814eSBarry Smith   SNESLineSearchGetReason - Gets the success/failure status of the last line search application
1633f40b411bSPeter Brune 
16342fe279fdSBarry Smith   Input Parameter:
1635*420bcc1bSBarry Smith . linesearch - the line search context
1636f40b411bSPeter Brune 
16372fe279fdSBarry Smith   Output Parameter:
1638422a814eSBarry Smith . result - The success or failure status
1639f40b411bSPeter Brune 
164020f4b53cSBarry Smith   Level: developer
164120f4b53cSBarry Smith 
1642f6dfbefdSBarry Smith   Note:
1643*420bcc1bSBarry Smith   This is typically called after `SNESLineSearchApply()` in order to determine if the line search failed
1644*420bcc1bSBarry Smith   (and set into the `SNES` convergence accordingly).
1645cd7522eaSPeter Brune 
1646*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetReason()`, `SNESLineSearchReason`
1647f40b411bSPeter Brune @*/
1648d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result)
1649d71ae5a4SJacob Faibussowitsch {
1650bf7f4e0aSPeter Brune   PetscFunctionBegin;
1651f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
16524f572ea9SToby Isaac   PetscAssertPointer(result, 2);
1653422a814eSBarry Smith   *result = linesearch->result;
16543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1655bf7f4e0aSPeter Brune }
1656bf7f4e0aSPeter Brune 
1657f40b411bSPeter Brune /*@
1658*420bcc1bSBarry Smith   SNESLineSearchSetReason - Sets the success/failure status of the line search application
1659f40b411bSPeter Brune 
1660f40b411bSPeter Brune   Input Parameters:
1661*420bcc1bSBarry Smith + linesearch - the line search context
1662422a814eSBarry Smith - result     - The success or failure status
1663f40b411bSPeter Brune 
166420f4b53cSBarry Smith   Level: developer
166520f4b53cSBarry Smith 
1666f6dfbefdSBarry Smith   Note:
1667*420bcc1bSBarry Smith   This is typically called in a `SNESLineSearchType` implementation of `SNESLineSearchApply()` or a `SNESLINESEARCHSHELL` implementation to set
1668cd7522eaSPeter Brune   the success or failure of the line search method.
1669cd7522eaSPeter Brune 
1670*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchReason`, `SNESLineSearchGetSResult()`
1671f40b411bSPeter Brune @*/
1672d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result)
1673d71ae5a4SJacob Faibussowitsch {
16746a388c36SPeter Brune   PetscFunctionBegin;
16755fd66863SKarl Rupp   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1676422a814eSBarry Smith   linesearch->result = result;
16773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16786a388c36SPeter Brune }
16796a388c36SPeter Brune 
1680ceaaa498SBarry Smith // PetscClangLinter pragma disable: -fdoc-param-list-func-parameter-documentation
16819bd66eb0SPeter Brune /*@C
1682f1c6b773SPeter Brune   SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.
16839bd66eb0SPeter Brune 
1684c3339decSBarry Smith   Logically Collective
1685f6dfbefdSBarry Smith 
16869bd66eb0SPeter Brune   Input Parameters:
1687ceaaa498SBarry Smith + linesearch  - the linesearch object
1688*420bcc1bSBarry Smith . projectfunc - function for projecting the function to the bounds, see `SNESLineSearchSetVIFunctions` for calling sequence
1689*420bcc1bSBarry Smith - normfunc    - function for computing the norm of an active set, see `SNESLineSearchSetVIFunctions ` for calling sequence
16909bd66eb0SPeter Brune 
1691f6dfbefdSBarry Smith   Level: advanced
16929bd66eb0SPeter Brune 
1693ceaaa498SBarry Smith   Notes:
169420f4b53cSBarry Smith   The VI solvers require projection of the solution to the feasible set.  `projectfunc` should implement this.
169520f4b53cSBarry Smith 
169620f4b53cSBarry Smith   The VI solvers require special evaluation of the function norm such that the norm is only calculated
169720f4b53cSBarry Smith   on the inactive set.  This should be implemented by `normfunc`.
169820f4b53cSBarry Smith 
1699*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetVIFunctions()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()`
17009bd66eb0SPeter Brune @*/
1701d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc)
1702d71ae5a4SJacob Faibussowitsch {
17039bd66eb0SPeter Brune   PetscFunctionBegin;
1704f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
17059bd66eb0SPeter Brune   if (projectfunc) linesearch->ops->viproject = projectfunc;
17069bd66eb0SPeter Brune   if (normfunc) linesearch->ops->vinorm = normfunc;
17073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17089bd66eb0SPeter Brune }
17099bd66eb0SPeter Brune 
17109bd66eb0SPeter Brune /*@C
1711f1c6b773SPeter Brune   SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.
17129bd66eb0SPeter Brune 
1713f6dfbefdSBarry Smith   Not Collective
1714f6dfbefdSBarry Smith 
1715f899ff85SJose E. Roman   Input Parameter:
1716f6dfbefdSBarry Smith . linesearch - the line search context, obtain with `SNESGetLineSearch()`
17179bd66eb0SPeter Brune 
17189bd66eb0SPeter Brune   Output Parameters:
17199bd66eb0SPeter Brune + projectfunc - function for projecting the function to the bounds
17209bd66eb0SPeter Brune - normfunc    - function for computing the norm of an active set
17219bd66eb0SPeter Brune 
1722f6dfbefdSBarry Smith   Level: advanced
17239bd66eb0SPeter Brune 
1724*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetVIFunctions()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()`
17259bd66eb0SPeter Brune @*/
1726d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc)
1727d71ae5a4SJacob Faibussowitsch {
17289bd66eb0SPeter Brune   PetscFunctionBegin;
17299bd66eb0SPeter Brune   if (projectfunc) *projectfunc = linesearch->ops->viproject;
17309bd66eb0SPeter Brune   if (normfunc) *normfunc = linesearch->ops->vinorm;
17313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17329bd66eb0SPeter Brune }
17339bd66eb0SPeter Brune 
1734bf7f4e0aSPeter Brune /*@C
1735*420bcc1bSBarry Smith   SNESLineSearchRegister - register a line search type `SNESLineSearchType`
1736ceaaa498SBarry Smith 
1737ceaaa498SBarry Smith   Input Parameters:
1738*420bcc1bSBarry Smith + sname    - name of the `SNESLineSearchType()`
1739ceaaa498SBarry Smith - function - the creation function for that type
1740ceaaa498SBarry Smith 
1741ceaaa498SBarry Smith   Calling sequence of `function`:
1742ceaaa498SBarry Smith . ls - the line search context
1743bf7f4e0aSPeter Brune 
1744bf7f4e0aSPeter Brune   Level: advanced
1745f6dfbefdSBarry Smith 
1746*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchSetType()`
1747bf7f4e0aSPeter Brune @*/
1748ceaaa498SBarry Smith PetscErrorCode SNESLineSearchRegister(const char sname[], PetscErrorCode (*function)(SNESLineSearch ls))
1749d71ae5a4SJacob Faibussowitsch {
1750bf7f4e0aSPeter Brune   PetscFunctionBegin;
17519566063dSJacob Faibussowitsch   PetscCall(SNESInitializePackage());
17529566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&SNESLineSearchList, sname, function));
17533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1754bf7f4e0aSPeter Brune }
1755