xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision 648c30bcb65f74c3cbd15d7a91a7ed7c1890e25b)
1af0996ceSBarry Smith #include <petsc/private/linesearchimpl.h> /*I "petscsnes.h" I*/
2bf7f4e0aSPeter Brune 
3f1c6b773SPeter Brune PetscBool         SNESLineSearchRegisterAllCalled = PETSC_FALSE;
40298fd71SBarry Smith PetscFunctionList SNESLineSearchList              = NULL;
5bf7f4e0aSPeter Brune 
6f1c6b773SPeter Brune PetscClassId  SNESLINESEARCH_CLASSID;
757a83faaSBarry Smith PetscLogEvent SNESLINESEARCH_Apply;
8bf7f4e0aSPeter Brune 
9dcf2fd19SBarry Smith /*@
10f6dfbefdSBarry Smith   SNESLineSearchMonitorCancel - Clears all the monitor functions for a `SNESLineSearch` object.
11dcf2fd19SBarry Smith 
12c3339decSBarry Smith   Logically Collective
13dcf2fd19SBarry Smith 
142fe279fdSBarry Smith   Input Parameter:
15f6dfbefdSBarry Smith . ls - the `SNESLineSearch` context
16dcf2fd19SBarry Smith 
17dcf2fd19SBarry Smith   Options Database Key:
18dcf2fd19SBarry Smith . -snes_linesearch_monitor_cancel - cancels all monitors that have been hardwired
19f6dfbefdSBarry Smith     into a code by calls to `SNESLineSearchMonitorSet()`, but does not cancel those
20dcf2fd19SBarry Smith     set via the options database
21dcf2fd19SBarry Smith 
2220f4b53cSBarry Smith   Level: advanced
2320f4b53cSBarry Smith 
24dcf2fd19SBarry Smith   Notes:
25f6dfbefdSBarry Smith   There is no way to clear one specific monitor from a `SNESLineSearch` object.
26dcf2fd19SBarry Smith 
27420bcc1bSBarry Smith   This does not clear the monitor set with `SNESLineSearchSetDefaultMonitor()` use `SNESLineSearchSetDefaultMonitor`(`ls`,`NULL`) to cancel it
28dcf2fd19SBarry Smith   that one.
29dcf2fd19SBarry Smith 
30420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorDefault()`, `SNESLineSearchMonitorSet()`
31dcf2fd19SBarry Smith @*/
32d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchMonitorCancel(SNESLineSearch ls)
33d71ae5a4SJacob Faibussowitsch {
34dcf2fd19SBarry Smith   PetscInt i;
35dcf2fd19SBarry Smith 
36dcf2fd19SBarry Smith   PetscFunctionBegin;
37dcf2fd19SBarry Smith   PetscValidHeaderSpecific(ls, SNESLINESEARCH_CLASSID, 1);
38dcf2fd19SBarry Smith   for (i = 0; i < ls->numbermonitors; i++) {
3948a46eb9SPierre Jolivet     if (ls->monitordestroy[i]) PetscCall((*ls->monitordestroy[i])(&ls->monitorcontext[i]));
40dcf2fd19SBarry Smith   }
41dcf2fd19SBarry Smith   ls->numbermonitors = 0;
423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
43dcf2fd19SBarry Smith }
44dcf2fd19SBarry Smith 
45dcf2fd19SBarry Smith /*@
46dcf2fd19SBarry Smith   SNESLineSearchMonitor - runs the user provided monitor routines, if they exist
47dcf2fd19SBarry Smith 
48c3339decSBarry Smith   Collective
49dcf2fd19SBarry Smith 
502fe279fdSBarry Smith   Input Parameter:
51dcf2fd19SBarry Smith . ls - the linesearch object
52dcf2fd19SBarry Smith 
532fe279fdSBarry Smith   Level: developer
542fe279fdSBarry Smith 
55f6dfbefdSBarry Smith   Note:
56420bcc1bSBarry Smith   This routine is called by the `SNESLineSearch` implementations.
57dcf2fd19SBarry Smith   It does not typically need to be called by the user.
58dcf2fd19SBarry Smith 
59420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorSet()`
60dcf2fd19SBarry Smith @*/
61d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchMonitor(SNESLineSearch ls)
62d71ae5a4SJacob Faibussowitsch {
63dcf2fd19SBarry Smith   PetscInt i, n = ls->numbermonitors;
64dcf2fd19SBarry Smith 
65dcf2fd19SBarry Smith   PetscFunctionBegin;
6648a46eb9SPierre Jolivet   for (i = 0; i < n; i++) PetscCall((*ls->monitorftns[i])(ls, ls->monitorcontext[i]));
673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
68dcf2fd19SBarry Smith }
69dcf2fd19SBarry Smith 
70dcf2fd19SBarry Smith /*@C
71dcf2fd19SBarry Smith   SNESLineSearchMonitorSet - Sets an ADDITIONAL function that is to be used at every
72dcf2fd19SBarry Smith   iteration of the nonlinear solver to display the iteration's
73dcf2fd19SBarry Smith   progress.
74dcf2fd19SBarry Smith 
75c3339decSBarry Smith   Logically Collective
76dcf2fd19SBarry Smith 
77dcf2fd19SBarry Smith   Input Parameters:
78f6dfbefdSBarry Smith + ls             - the `SNESLineSearch` context
79dcf2fd19SBarry Smith . f              - the monitor function
80420bcc1bSBarry Smith . mctx           - [optional] user-defined context for private data for the monitor routine (use `NULL` if no context is desired)
81420bcc1bSBarry Smith - monitordestroy - [optional] routine that frees monitor context (may be `NULL`)
82420bcc1bSBarry Smith 
83420bcc1bSBarry Smith   Calling sequence of `f`:
84420bcc1bSBarry Smith + ls   - the `SNESLineSearch` context
85420bcc1bSBarry Smith - mctx - [optional] user-defined context for private data for the monitor routine
86420bcc1bSBarry Smith 
87420bcc1bSBarry Smith   Calling sequence of `monitordestroy`:
88420bcc1bSBarry 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 
97420bcc1bSBarry Smith   Fortran Note:
98f6dfbefdSBarry Smith   Only a single monitor function can be set for each `SNESLineSearch` object
99dcf2fd19SBarry Smith 
100420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorDefault()`, `SNESLineSearchMonitorCancel()`
101dcf2fd19SBarry Smith @*/
102420bcc1bSBarry 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:
126420bcc1bSBarry Smith + ls - the `SNESLineSearch` object
127f6dfbefdSBarry Smith - vf - the context for the monitor, in this case it is an `PetscViewerAndFormat`
128dcf2fd19SBarry Smith 
129420bcc1bSBarry Smith   Options Database Key:
130420bcc1bSBarry Smith . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine
131dcf2fd19SBarry Smith 
132420bcc1bSBarry Smith   Level: developer
133420bcc1bSBarry Smith 
134420bcc1bSBarry Smith   This is not normally called directly but is passed to `SNESLineSearchMonitorSet()`
135420bcc1bSBarry Smith 
136420bcc1bSBarry 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 /*@
157420bcc1bSBarry 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:
170420bcc1bSBarry Smith   The preferred calling sequence is to use `SNESGetLineSearch()` to acquire the `SNESLineSearch` instance
171f6dfbefdSBarry Smith   already associated with the `SNES`.
172f40b411bSPeter Brune 
173420bcc1bSBarry 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());
182f5af7f23SKarl Rupp 
1839566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(linesearch, SNESLINESEARCH_CLASSID, "SNESLineSearch", "Linesearch", "SNESLineSearch", comm, SNESLineSearchDestroy, SNESLineSearchView));
1840298fd71SBarry Smith   linesearch->vec_sol_new  = NULL;
1850298fd71SBarry Smith   linesearch->vec_func_new = NULL;
1860298fd71SBarry Smith   linesearch->vec_sol      = NULL;
1870298fd71SBarry Smith   linesearch->vec_func     = NULL;
1880298fd71SBarry Smith   linesearch->vec_update   = NULL;
1899bd66eb0SPeter Brune 
190bf7f4e0aSPeter Brune   linesearch->lambda       = 1.0;
191bf7f4e0aSPeter Brune   linesearch->fnorm        = 1.0;
192bf7f4e0aSPeter Brune   linesearch->ynorm        = 1.0;
193bf7f4e0aSPeter Brune   linesearch->xnorm        = 1.0;
194422a814eSBarry Smith   linesearch->result       = SNES_LINESEARCH_SUCCEEDED;
195bf7f4e0aSPeter Brune   linesearch->norms        = PETSC_TRUE;
196bf7f4e0aSPeter Brune   linesearch->keeplambda   = PETSC_FALSE;
197bf7f4e0aSPeter Brune   linesearch->damping      = 1.0;
198bf7f4e0aSPeter Brune   linesearch->maxstep      = 1e8;
199bf7f4e0aSPeter Brune   linesearch->steptol      = 1e-12;
200516fe3c3SPeter Brune   linesearch->rtol         = 1e-8;
201516fe3c3SPeter Brune   linesearch->atol         = 1e-15;
202516fe3c3SPeter Brune   linesearch->ltol         = 1e-8;
2030298fd71SBarry Smith   linesearch->precheckctx  = NULL;
2040298fd71SBarry Smith   linesearch->postcheckctx = NULL;
20559405d5eSPeter Brune   linesearch->max_its      = 1;
206bf7f4e0aSPeter Brune   linesearch->setupcalled  = PETSC_FALSE;
2073add74b1SBarry Smith   linesearch->monitor      = NULL;
208bf7f4e0aSPeter Brune   *outlinesearch           = linesearch;
2093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
210bf7f4e0aSPeter Brune }
211bf7f4e0aSPeter Brune 
212f40b411bSPeter Brune /*@
21378bcb3b5SPeter Brune   SNESLineSearchSetUp - Prepares the line search for being applied by allocating
21478bcb3b5SPeter Brune   any required vectors.
215f40b411bSPeter Brune 
216c3339decSBarry Smith   Collective
217f40b411bSPeter Brune 
2182fe279fdSBarry Smith   Input Parameter:
219f6dfbefdSBarry Smith . linesearch - The `SNESLineSearch` instance.
220f40b411bSPeter Brune 
22120f4b53cSBarry Smith   Level: advanced
22220f4b53cSBarry Smith 
223f6dfbefdSBarry Smith   Note:
224f6dfbefdSBarry Smith   For most cases, this needn't be called by users or outside of `SNESLineSearchApply()`.
225cd7522eaSPeter Brune   The only current case where this is called outside of this is for the VI
22678bcb3b5SPeter Brune   solvers, which modify the solution and work vectors before the first call
227f6dfbefdSBarry Smith   of `SNESLineSearchApply()`, requiring the `SNESLineSearch` work vectors to be
228cd7522eaSPeter Brune   allocated upfront.
229cd7522eaSPeter Brune 
230420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchReset()`
231f40b411bSPeter Brune @*/
232d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch)
233d71ae5a4SJacob Faibussowitsch {
234bf7f4e0aSPeter Brune   PetscFunctionBegin;
23548a46eb9SPierre Jolivet   if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC));
236bf7f4e0aSPeter Brune   if (!linesearch->setupcalled) {
23748a46eb9SPierre Jolivet     if (!linesearch->vec_sol_new) PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new));
23848a46eb9SPierre Jolivet     if (!linesearch->vec_func_new) PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new));
239213acdd3SPierre Jolivet     PetscTryTypeMethod(linesearch, setup);
2409566063dSJacob Faibussowitsch     if (!linesearch->ops->snesfunc) PetscCall(SNESLineSearchSetFunction(linesearch, SNESComputeFunction));
241bf7f4e0aSPeter Brune     linesearch->lambda      = linesearch->damping;
242bf7f4e0aSPeter Brune     linesearch->setupcalled = PETSC_TRUE;
243bf7f4e0aSPeter Brune   }
2443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
245bf7f4e0aSPeter Brune }
246bf7f4e0aSPeter Brune 
247f40b411bSPeter Brune /*@
248f6dfbefdSBarry Smith   SNESLineSearchReset - Undoes the `SNESLineSearchSetUp()` and deletes any `Vec`s or `Mat`s allocated by the line search.
249f40b411bSPeter Brune 
250c3339decSBarry Smith   Collective
251f40b411bSPeter Brune 
2522fe279fdSBarry Smith   Input Parameter:
253f6dfbefdSBarry Smith . linesearch - The `SNESLineSearch` instance.
254f40b411bSPeter Brune 
25520f4b53cSBarry Smith   Level: developer
25620f4b53cSBarry Smith 
257f6dfbefdSBarry Smith   Note:
258f6dfbefdSBarry Smith   Usually only called by `SNESReset()`
259f190f2fcSBarry Smith 
260420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchSetUp()`
261f40b411bSPeter Brune @*/
262d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch)
263d71ae5a4SJacob Faibussowitsch {
264bf7f4e0aSPeter Brune   PetscFunctionBegin;
265213acdd3SPierre Jolivet   PetscTryTypeMethod(linesearch, reset);
266f5af7f23SKarl Rupp 
2679566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&linesearch->vec_sol_new));
2689566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&linesearch->vec_func_new));
269bf7f4e0aSPeter Brune 
2709566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(linesearch->nwork, &linesearch->work));
271f5af7f23SKarl Rupp 
272bf7f4e0aSPeter Brune   linesearch->nwork       = 0;
273bf7f4e0aSPeter Brune   linesearch->setupcalled = PETSC_FALSE;
2743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
275bf7f4e0aSPeter Brune }
276bf7f4e0aSPeter Brune 
277ed07d7d7SPeter Brune /*@C
278f6dfbefdSBarry Smith   SNESLineSearchSetFunction - Sets the function evaluation used by the `SNES` line search
279f6dfbefdSBarry Smith   `
280e4094ef1SJacob Faibussowitsch 
281ed07d7d7SPeter Brune   Input Parameters:
282e4094ef1SJacob Faibussowitsch + linesearch - the `SNESLineSearch` context
283e4094ef1SJacob Faibussowitsch - func       - function evaluation routine, this is usually the function provided with `SNESSetFunction()`
284ed07d7d7SPeter Brune 
285420bcc1bSBarry Smith   Calling sequence of `func`:
286420bcc1bSBarry Smith + snes - the `SNES` with which the `SNESLineSearch` context is associated with
287420bcc1bSBarry Smith . x    - the input vector
288420bcc1bSBarry Smith - f    - the computed value of the function
289420bcc1bSBarry Smith 
290ed07d7d7SPeter Brune   Level: developer
291ed07d7d7SPeter Brune 
292420bcc1bSBarry Smith   Note:
293420bcc1bSBarry Smith   By default the `SNESLineSearch` uses the function provided by `SNESSetFunction()` so this is rarely needed
294420bcc1bSBarry Smith 
295420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESSetFunction()`
296ed07d7d7SPeter Brune @*/
297420bcc1bSBarry Smith PetscErrorCode SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES snes, Vec x, Vec f))
298d71ae5a4SJacob Faibussowitsch {
299ed07d7d7SPeter Brune   PetscFunctionBegin;
300ed07d7d7SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
301ed07d7d7SPeter Brune   linesearch->ops->snesfunc = func;
3023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
303ed07d7d7SPeter Brune }
304ed07d7d7SPeter Brune 
30586d74e61SPeter Brune /*@C
306420bcc1bSBarry Smith   SNESLineSearchSetPreCheck - Sets a function that is called after the initial search direction has been computed but
307420bcc1bSBarry Smith   before the line search routine has been applied. Allows adjusting the result of (usually a linear solve) that
308f190f2fcSBarry Smith   determined the search direction.
30986d74e61SPeter Brune 
310c3339decSBarry Smith   Logically Collective
31186d74e61SPeter Brune 
31286d74e61SPeter Brune   Input Parameters:
313f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context
314420bcc1bSBarry Smith . func       - [optional] function evaluation routine
31520f4b53cSBarry Smith - ctx        - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
31686d74e61SPeter Brune 
317420bcc1bSBarry Smith   Calling sequence of `func`:
318420bcc1bSBarry Smith + ls        - the `SNESLineSearch` context
319420bcc1bSBarry Smith . x         - the current solution
320a678f235SPierre Jolivet . d         - the current search direction
321420bcc1bSBarry Smith . changed_d - indicates if the search direction has been changed
322420bcc1bSBarry Smith - ctx       - the context passed to `SNESLineSearchSetPreCheck()`
323420bcc1bSBarry Smith 
32486d74e61SPeter Brune   Level: intermediate
32586d74e61SPeter Brune 
326f6dfbefdSBarry Smith   Note:
327420bcc1bSBarry Smith   Use `SNESLineSearchSetPostCheck()` to change the step after the line search is complete.
328f0b84518SBarry Smith 
329f0b84518SBarry Smith   Use `SNESVISetVariableBounds()` and `SNESVISetComputeVariableBounds()` to cause `SNES` to automatically control the ranges of variables allowed.
330f0b84518SBarry Smith 
331420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()`,
332869ba2dcSBarry Smith           `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`
333f0b84518SBarry Smith 
33486d74e61SPeter Brune @*/
335420bcc1bSBarry Smith PetscErrorCode SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch ls, Vec x, Vec d, PetscBool *changed_d, void *ctx), void *ctx)
336d71ae5a4SJacob Faibussowitsch {
3379bd66eb0SPeter Brune   PetscFunctionBegin;
338f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
339f190f2fcSBarry Smith   if (func) linesearch->ops->precheck = func;
34086d74e61SPeter Brune   if (ctx) linesearch->precheckctx = ctx;
3413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
34286d74e61SPeter Brune }
34386d74e61SPeter Brune 
34486d74e61SPeter Brune /*@C
34553deb899SBarry Smith   SNESLineSearchGetPreCheck - Gets the pre-check function for the line search routine.
34686d74e61SPeter Brune 
347f899ff85SJose E. Roman   Input Parameter:
348f6dfbefdSBarry Smith . linesearch - the `SNESLineSearch` context
34986d74e61SPeter Brune 
35086d74e61SPeter Brune   Output Parameters:
351420bcc1bSBarry Smith + func - [optional] function evaluation routine,  for calling sequence see `SNESLineSearchSetPreCheck()`
35220f4b53cSBarry Smith - ctx  - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
35386d74e61SPeter Brune 
35486d74e61SPeter Brune   Level: intermediate
35586d74e61SPeter Brune 
356420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchSetPostCheck()`
35786d74e61SPeter Brune @*/
358d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch, Vec, Vec, PetscBool *, void *), void **ctx)
359d71ae5a4SJacob Faibussowitsch {
3609bd66eb0SPeter Brune   PetscFunctionBegin;
361f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
362f190f2fcSBarry Smith   if (func) *func = linesearch->ops->precheck;
36386d74e61SPeter Brune   if (ctx) *ctx = linesearch->precheckctx;
3643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
36586d74e61SPeter Brune }
36686d74e61SPeter Brune 
36786d74e61SPeter Brune /*@C
368f190f2fcSBarry Smith   SNESLineSearchSetPostCheck - Sets a user function that is called after the line search has been applied to determine the step
369f190f2fcSBarry Smith   direction and length. Allows the user a chance to change or override the decision of the line search routine
37086d74e61SPeter Brune 
37120f4b53cSBarry Smith   Logically Collective
37286d74e61SPeter Brune 
37386d74e61SPeter Brune   Input Parameters:
374f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context
375420bcc1bSBarry Smith . func       - [optional] function evaluation routine
37620f4b53cSBarry Smith - ctx        - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
37786d74e61SPeter Brune 
378420bcc1bSBarry Smith   Calling sequence of `func`:
379420bcc1bSBarry Smith + ls        - the `SNESLineSearch` context
380420bcc1bSBarry Smith . x         - the current solution
381a678f235SPierre Jolivet . d         - the current search direction
382a678f235SPierre Jolivet . w         - $ w = x + lambda*d $ for some lambda
383420bcc1bSBarry Smith . changed_d - indicates if the search direction `d` has been changed
384420bcc1bSBarry Smith . changed_w - indicates `w` has been changed
385420bcc1bSBarry Smith - ctx       - the context passed to `SNESLineSearchSetPreCheck()`
386420bcc1bSBarry Smith 
38786d74e61SPeter Brune   Level: intermediate
38886d74e61SPeter Brune 
389f0b84518SBarry Smith   Notes:
3906b095a85SStefano Zampini   Use `SNESLineSearchSetPreCheck()` to change the step before the line search is completed.
3916b095a85SStefano Zampini   The calling sequence of the callback does not contain the current scaling factor. To access the value, use `SNESLineSearchGetLambda()`.
392f0b84518SBarry Smith 
393f0b84518SBarry Smith   Use `SNESVISetVariableBounds()` and `SNESVISetComputeVariableBounds()` to cause `SNES` to automatically control the ranges of variables allowed.
394f0b84518SBarry Smith 
395420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchGetPostCheck()`,
396e4094ef1SJacob Faibussowitsch           `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`
39786d74e61SPeter Brune @*/
398420bcc1bSBarry 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)
399d71ae5a4SJacob Faibussowitsch {
40086d74e61SPeter Brune   PetscFunctionBegin;
401f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
402f190f2fcSBarry Smith   if (func) linesearch->ops->postcheck = func;
40386d74e61SPeter Brune   if (ctx) linesearch->postcheckctx = ctx;
4043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
40586d74e61SPeter Brune }
40686d74e61SPeter Brune 
40786d74e61SPeter Brune /*@C
408f1c6b773SPeter Brune   SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine.
40986d74e61SPeter Brune 
410f899ff85SJose E. Roman   Input Parameter:
411f6dfbefdSBarry Smith . linesearch - the `SNESLineSearch` context
41286d74e61SPeter Brune 
41386d74e61SPeter Brune   Output Parameters:
414420bcc1bSBarry Smith + func - [optional] function evaluation routine, see for the calling sequence `SNESLineSearchSetPostCheck()`
41520f4b53cSBarry Smith - ctx  - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
41686d74e61SPeter Brune 
41786d74e61SPeter Brune   Level: intermediate
41886d74e61SPeter Brune 
419420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`
42086d74e61SPeter Brune @*/
421d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx)
422d71ae5a4SJacob Faibussowitsch {
4239bd66eb0SPeter Brune   PetscFunctionBegin;
424f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
425f190f2fcSBarry Smith   if (func) *func = linesearch->ops->postcheck;
42686d74e61SPeter Brune   if (ctx) *ctx = linesearch->postcheckctx;
4273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
42886d74e61SPeter Brune }
42986d74e61SPeter Brune 
430f40b411bSPeter Brune /*@
431f1c6b773SPeter Brune   SNESLineSearchPreCheck - Prepares the line search for being applied.
432f40b411bSPeter Brune 
433c3339decSBarry Smith   Logically Collective
434f40b411bSPeter Brune 
435f40b411bSPeter Brune   Input Parameters:
4367b1df9c1SPeter Brune + linesearch - The linesearch instance.
4377b1df9c1SPeter Brune . X          - The current solution
4387b1df9c1SPeter Brune - Y          - The step direction
439f40b411bSPeter Brune 
4402fe279fdSBarry Smith   Output Parameter:
441420bcc1bSBarry Smith . changed - Indicator that the precheck routine has changed `Y`
442f40b411bSPeter Brune 
443f0b84518SBarry Smith   Level: advanced
444f40b411bSPeter Brune 
445f6dfbefdSBarry Smith   Note:
446420bcc1bSBarry Smith   This calls any function provided with `SNESLineSearchSetPreCheck()` and is called automatically inside the line search routines
447f6dfbefdSBarry Smith 
448420bcc1bSBarry Smith   Developer Note:
449420bcc1bSBarry Smith   The use of `PetscObjectGetState()` would eliminate the need for the `changed` argument to be provided
450420bcc1bSBarry Smith 
451420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`,
452fe8e7dddSPierre Jolivet           `SNESLineSearchGetPostCheck()`
453f40b411bSPeter Brune @*/
454d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed)
455d71ae5a4SJacob Faibussowitsch {
456bf7f4e0aSPeter Brune   PetscFunctionBegin;
457bf7f4e0aSPeter Brune   *changed = PETSC_FALSE;
4586b2b7091SBarry Smith   if (linesearch->ops->precheck) {
459dbbe0bcdSBarry Smith     PetscUseTypeMethod(linesearch, precheck, X, Y, changed, linesearch->precheckctx);
46038bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch, *changed, 4);
461bf7f4e0aSPeter Brune   }
4623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
463bf7f4e0aSPeter Brune }
464bf7f4e0aSPeter Brune 
465f40b411bSPeter Brune /*@
466ef46b1a6SStefano Zampini   SNESLineSearchPostCheck - Hook to modify step direction or updated solution after a successful linesearch
467f40b411bSPeter Brune 
468c3339decSBarry Smith   Logically Collective
469f40b411bSPeter Brune 
470f40b411bSPeter Brune   Input Parameters:
4717b1df9c1SPeter Brune + linesearch - The line search context
4727b1df9c1SPeter Brune . X          - The last solution
4737b1df9c1SPeter Brune . Y          - The step direction
4746b095a85SStefano Zampini - W          - The updated solution, `W = X - lambda * Y` for some lambda
475f40b411bSPeter Brune 
476f40b411bSPeter Brune   Output Parameters:
4776b095a85SStefano Zampini + changed_Y - Indicator if the direction `Y` has been changed.
478420bcc1bSBarry Smith - changed_W - Indicator if the new candidate solution `W` has been changed.
479f40b411bSPeter Brune 
48020f4b53cSBarry Smith   Level: developer
48120f4b53cSBarry Smith 
482f6dfbefdSBarry Smith   Note:
483420bcc1bSBarry Smith   This calls any function provided with `SNESLineSearchSetPostCheck()` and is called automatically inside the line search routines
484f6dfbefdSBarry Smith 
485420bcc1bSBarry Smith   Developer Note:
486420bcc1bSBarry Smith   The use of `PetscObjectGetState()` would eliminate the need for the `changed_Y` and `changed_W` arguments to be provided
487420bcc1bSBarry Smith 
488420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPrecheck()`, `SNESLineSearchGetPrecheck()`
489f40b411bSPeter Brune @*/
490d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W)
491d71ae5a4SJacob Faibussowitsch {
492bf7f4e0aSPeter Brune   PetscFunctionBegin;
493bf7f4e0aSPeter Brune   *changed_Y = PETSC_FALSE;
494bf7f4e0aSPeter Brune   *changed_W = PETSC_FALSE;
4956b2b7091SBarry Smith   if (linesearch->ops->postcheck) {
496dbbe0bcdSBarry Smith     PetscUseTypeMethod(linesearch, postcheck, X, Y, W, changed_Y, changed_W, linesearch->postcheckctx);
49738bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch, *changed_Y, 5);
49838bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch, *changed_W, 6);
49986d74e61SPeter Brune   }
5003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
50186d74e61SPeter Brune }
50286d74e61SPeter Brune 
50386d74e61SPeter Brune /*@C
5041d27aa22SBarry Smith   SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration {cite}`hindmarsh1996time`
50586d74e61SPeter Brune 
506c3339decSBarry Smith   Logically Collective
50786d74e61SPeter Brune 
5084165533cSJose E. Roman   Input Parameters:
509420bcc1bSBarry Smith + linesearch - the line search context
51086d74e61SPeter Brune . X          - base state for this step
511907376e6SBarry Smith - ctx        - context for this function
51286d74e61SPeter Brune 
51397bb3fdcSJose E. Roman   Input/Output Parameter:
51497bb3fdcSJose E. Roman . Y - correction, possibly modified
51597bb3fdcSJose E. Roman 
51697bb3fdcSJose E. Roman   Output Parameter:
517420bcc1bSBarry Smith . changed - flag indicating that `Y` was modified
51886d74e61SPeter Brune 
519420bcc1bSBarry Smith   Options Database Keys:
520cd7522eaSPeter Brune + -snes_linesearch_precheck_picard       - activate this routine
521cd7522eaSPeter Brune - -snes_linesearch_precheck_picard_angle - angle
52286d74e61SPeter Brune 
52386d74e61SPeter Brune   Level: advanced
52486d74e61SPeter Brune 
52586d74e61SPeter Brune   Notes:
526f6dfbefdSBarry Smith   This function should be passed to `SNESLineSearchSetPreCheck()`
52786d74e61SPeter Brune 
52886d74e61SPeter Brune   The justification for this method involves the linear convergence of a Picard iteration
5291d27aa22SBarry Smith   so the Picard linearization should be provided in place of the "Jacobian"  {cite}`hindmarsh1996time`. This correction
53086d74e61SPeter Brune   is generally not useful when using a Newton linearization.
53186d74e61SPeter Brune 
532420bcc1bSBarry Smith   Developer Note:
533420bcc1bSBarry Smith   The use of `PetscObjectGetState()` would eliminate the need for the `changed` argument to be provided
534420bcc1bSBarry Smith 
535420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESSetPicard()`, `SNESGetLineSearch()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchSetPostCheck()`
53686d74e61SPeter Brune @*/
537d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed, void *ctx)
538d71ae5a4SJacob Faibussowitsch {
53986d74e61SPeter Brune   PetscReal   angle = *(PetscReal *)linesearch->precheckctx;
54086d74e61SPeter Brune   Vec         Ylast;
54186d74e61SPeter Brune   PetscScalar dot;
54286d74e61SPeter Brune   PetscInt    iter;
54386d74e61SPeter Brune   PetscReal   ynorm, ylastnorm, theta, angle_radians;
54486d74e61SPeter Brune   SNES        snes;
54586d74e61SPeter Brune 
54686d74e61SPeter Brune   PetscFunctionBegin;
5479566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
5489566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)snes, "SNESLineSearchPreCheckPicard_Ylast", (PetscObject *)&Ylast));
54986d74e61SPeter Brune   if (!Ylast) {
5509566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(Y, &Ylast));
5519566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)snes, "SNESLineSearchPreCheckPicard_Ylast", (PetscObject)Ylast));
5529566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)Ylast));
55386d74e61SPeter Brune   }
5549566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &iter));
55586d74e61SPeter Brune   if (iter < 2) {
5569566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y, Ylast));
55786d74e61SPeter Brune     *changed = PETSC_FALSE;
5583ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
55986d74e61SPeter Brune   }
56086d74e61SPeter Brune 
5619566063dSJacob Faibussowitsch   PetscCall(VecDot(Y, Ylast, &dot));
5629566063dSJacob Faibussowitsch   PetscCall(VecNorm(Y, NORM_2, &ynorm));
5639566063dSJacob Faibussowitsch   PetscCall(VecNorm(Ylast, NORM_2, &ylastnorm));
56486d74e61SPeter Brune   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
565255453a1SBarry Smith   theta         = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm), -1.0, 1.0));
56686d74e61SPeter Brune   angle_radians = angle * PETSC_PI / 180.;
56786d74e61SPeter Brune   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
56886d74e61SPeter Brune     /* Modify the step Y */
56986d74e61SPeter Brune     PetscReal alpha, ydiffnorm;
5709566063dSJacob Faibussowitsch     PetscCall(VecAXPY(Ylast, -1.0, Y));
5719566063dSJacob Faibussowitsch     PetscCall(VecNorm(Ylast, NORM_2, &ydiffnorm));
572f85e2ce2SBarry Smith     alpha = (ydiffnorm > .001 * ylastnorm) ? ylastnorm / ydiffnorm : 1000.0;
5739566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y, Ylast));
5749566063dSJacob Faibussowitsch     PetscCall(VecScale(Y, alpha));
5759566063dSJacob 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));
576a47ec85fSBarry Smith     *changed = PETSC_TRUE;
57786d74e61SPeter Brune   } else {
5789566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Angle %14.12e degrees exceeds threshold %14.12e, no correction applied\n", (double)(theta * 180. / PETSC_PI), (double)angle));
5799566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y, Ylast));
58086d74e61SPeter Brune     *changed = PETSC_FALSE;
581bf7f4e0aSPeter Brune   }
5823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
583bf7f4e0aSPeter Brune }
584bf7f4e0aSPeter Brune 
585f40b411bSPeter Brune /*@
586cd7522eaSPeter Brune   SNESLineSearchApply - Computes the line-search update.
587f40b411bSPeter Brune 
588c3339decSBarry Smith   Collective
589f40b411bSPeter Brune 
5906b095a85SStefano Zampini   Input Parameter:
5916b095a85SStefano Zampini . linesearch - The line search context
592f40b411bSPeter Brune 
5936b867d5aSJose E. Roman   Input/Output Parameters:
5946b867d5aSJose E. Roman + X     - The current solution, on output the new solution
595420bcc1bSBarry Smith . F     - The current function value, on output the new function value at the solution value `X`
5961bd63e3eSSatish Balay . fnorm - The current norm of `F`, on output the new norm of `F`
5976b095a85SStefano Zampini - Y     - The current search direction, on output the direction determined by the linesearch, i.e. Xnew = Xold - lambda*Y
598f40b411bSPeter Brune 
599cd7522eaSPeter Brune   Options Database Keys:
6000b00b554SBarry Smith + -snes_linesearch_type                - basic (or equivalently none), bt, l2, cp, nleqerr, shell
601dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches
6021fe24845SBarry Smith . -snes_linesearch_damping             - The linesearch damping parameter, default is 1.0 (no damping)
6031fe24845SBarry Smith . -snes_linesearch_norms               - Turn on/off the linesearch norms computation (SNESLineSearchSetComputeNorms())
6043c7d6663SPeter Brune . -snes_linesearch_keeplambda          - Keep the previous search length as the initial guess
6053c7d6663SPeter Brune - -snes_linesearch_max_it              - The number of iterations for iterative line searches
606cd7522eaSPeter Brune 
607e4094ef1SJacob Faibussowitsch   Level: intermediate
60820f4b53cSBarry Smith 
609cd7522eaSPeter Brune   Notes:
610f6dfbefdSBarry Smith   This is typically called from within a `SNESSolve()` implementation in order to
611f6dfbefdSBarry Smith   help with convergence of the nonlinear method.  Various `SNES` types use line searches
612cd7522eaSPeter Brune   in different ways, but the overarching theme is that a line search is used to determine
613cd7522eaSPeter Brune   an optimal damping parameter of a step at each iteration of the method.  Each
614f6dfbefdSBarry Smith   application of the line search may invoke `SNESComputeFunction()` several times, and
615cd7522eaSPeter Brune   therefore may be fairly expensive.
616cd7522eaSPeter Brune 
6171bd63e3eSSatish Balay .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchGetLambda()`, `SNESLineSearchPreCheck()`, `SNESLineSearchPostCheck()`, `SNESSolve()`, `SNESComputeFunction()`, `SNESLineSearchSetComputeNorms()`,
618db781477SPatrick Sanan           `SNESLineSearchType`, `SNESLineSearchSetType()`
619f40b411bSPeter Brune @*/
620d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal *fnorm, Vec Y)
621d71ae5a4SJacob Faibussowitsch {
622bf388a1fSBarry Smith   PetscFunctionBegin;
623f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
624bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
625bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
626064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(Y, VEC_CLASSID, 5);
627bf7f4e0aSPeter Brune 
628422a814eSBarry Smith   linesearch->result = SNES_LINESEARCH_SUCCEEDED;
629bf7f4e0aSPeter Brune 
630bf7f4e0aSPeter Brune   linesearch->vec_sol    = X;
631bf7f4e0aSPeter Brune   linesearch->vec_update = Y;
632bf7f4e0aSPeter Brune   linesearch->vec_func   = F;
633bf7f4e0aSPeter Brune 
6349566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchSetUp(linesearch));
635bf7f4e0aSPeter Brune 
636f5af7f23SKarl Rupp   if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */
637bf7f4e0aSPeter Brune 
638f5af7f23SKarl Rupp   if (fnorm) linesearch->fnorm = *fnorm;
63948a46eb9SPierre Jolivet   else PetscCall(VecNorm(F, NORM_2, &linesearch->fnorm));
640bf7f4e0aSPeter Brune 
6419566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(SNESLINESEARCH_Apply, linesearch, X, F, Y));
642bf7f4e0aSPeter Brune 
643dbbe0bcdSBarry Smith   PetscUseTypeMethod(linesearch, apply);
644bf7f4e0aSPeter Brune 
6459566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(SNESLINESEARCH_Apply, linesearch, X, F, Y));
646bf7f4e0aSPeter Brune 
647f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
6483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
649bf7f4e0aSPeter Brune }
650bf7f4e0aSPeter Brune 
651f40b411bSPeter Brune /*@
652f1c6b773SPeter Brune   SNESLineSearchDestroy - Destroys the line search instance.
653f40b411bSPeter Brune 
654c3339decSBarry Smith   Collective
655f40b411bSPeter Brune 
6562fe279fdSBarry Smith   Input Parameter:
6578e557f58SPeter Brune . linesearch - The line search context
658f40b411bSPeter Brune 
65984238204SBarry Smith   Level: developer
660f40b411bSPeter Brune 
661420bcc1bSBarry Smith   Note:
662420bcc1bSBarry Smith   The line search in `SNES` is automatically called on `SNESDestroy()` so this call is rarely needed
663420bcc1bSBarry Smith 
664420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchReset()`, `SNESDestroy()`
665f40b411bSPeter Brune @*/
666d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchDestroy(SNESLineSearch *linesearch)
667d71ae5a4SJacob Faibussowitsch {
668bf7f4e0aSPeter Brune   PetscFunctionBegin;
6693ba16761SJacob Faibussowitsch   if (!*linesearch) PetscFunctionReturn(PETSC_SUCCESS);
670f4f49eeaSPierre Jolivet   PetscValidHeaderSpecific(*linesearch, SNESLINESEARCH_CLASSID, 1);
671f4f49eeaSPierre Jolivet   if (--((PetscObject)*linesearch)->refct > 0) {
6729371c9d4SSatish Balay     *linesearch = NULL;
6733ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
6749371c9d4SSatish Balay   }
6759566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsViewOff((PetscObject)*linesearch));
6769566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchReset(*linesearch));
6773ba16761SJacob Faibussowitsch   PetscTryTypeMethod(*linesearch, destroy);
678*648c30bcSBarry Smith   PetscCall(PetscViewerDestroy(&(*linesearch)->monitor));
679f4f49eeaSPierre Jolivet   PetscCall(SNESLineSearchMonitorCancel(*linesearch));
6809566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(linesearch));
6813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
682bf7f4e0aSPeter Brune }
683bf7f4e0aSPeter Brune 
684f40b411bSPeter Brune /*@
685dcf2fd19SBarry Smith   SNESLineSearchSetDefaultMonitor - Turns on/off printing useful information and debugging output about the line search.
686bf7f4e0aSPeter Brune 
687c3339decSBarry Smith   Logically Collective
688f6dfbefdSBarry Smith 
689bf7f4e0aSPeter Brune   Input Parameters:
690dcf2fd19SBarry Smith + linesearch - the linesearch object
69120f4b53cSBarry Smith - viewer     - an `PETSCVIEWERASCII` `PetscViewer` or `NULL` to turn off monitor
692bf7f4e0aSPeter Brune 
693f6dfbefdSBarry Smith   Options Database Key:
694dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - enables the monitor
695bf7f4e0aSPeter Brune 
696bf7f4e0aSPeter Brune   Level: intermediate
697bf7f4e0aSPeter Brune 
698e4094ef1SJacob Faibussowitsch   Developer Notes:
699f6dfbefdSBarry Smith   This monitor is implemented differently than the other line search monitors that are set with
700f6dfbefdSBarry Smith   `SNESLineSearchMonitorSet()` since it is called in many locations of the line search routines to display aspects of the
701d12e167eSBarry Smith   line search that are not visible to the other monitors.
702dcf2fd19SBarry Smith 
703420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `PETSCVIEWERASCII`, `SNESGetLineSearch()`, `SNESLineSearchGetDefaultMonitor()`, `PetscViewer`, `SNESLineSearchSetMonitor()`,
704f6dfbefdSBarry Smith           `SNESLineSearchMonitorSetFromOptions()`
705bf7f4e0aSPeter Brune @*/
706d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetDefaultMonitor(SNESLineSearch linesearch, PetscViewer viewer)
707d71ae5a4SJacob Faibussowitsch {
708bf7f4e0aSPeter Brune   PetscFunctionBegin;
709*648c30bcSBarry Smith   PetscCall(PetscViewerDestroy(&linesearch->monitor));
710dcf2fd19SBarry Smith   linesearch->monitor = viewer;
7113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
712bf7f4e0aSPeter Brune }
713bf7f4e0aSPeter Brune 
714f40b411bSPeter Brune /*@
715420bcc1bSBarry Smith   SNESLineSearchGetDefaultMonitor - Gets the `PetscViewer` instance for the default line search monitor that is turned on with `SNESLineSearchSetDefaultMonitor()`
716f6dfbefdSBarry Smith 
717c3339decSBarry Smith   Logically Collective
7186a388c36SPeter Brune 
719f190f2fcSBarry Smith   Input Parameter:
720420bcc1bSBarry Smith . linesearch - the line search context
721f40b411bSPeter Brune 
722f190f2fcSBarry Smith   Output Parameter:
7238e557f58SPeter Brune . monitor - monitor context
724f40b411bSPeter Brune 
725f40b411bSPeter Brune   Level: intermediate
726f40b411bSPeter Brune 
727420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchSetDefaultMonitor()`, `PetscViewer`
728f40b411bSPeter Brune @*/
729d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetDefaultMonitor(SNESLineSearch linesearch, PetscViewer *monitor)
730d71ae5a4SJacob Faibussowitsch {
7316a388c36SPeter Brune   PetscFunctionBegin;
732f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
7336a388c36SPeter Brune   *monitor = linesearch->monitor;
7343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7356a388c36SPeter Brune }
7366a388c36SPeter Brune 
737dcf2fd19SBarry Smith /*@C
738420bcc1bSBarry Smith   SNESLineSearchMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated in the options database
739dcf2fd19SBarry Smith 
740c3339decSBarry Smith   Collective
741dcf2fd19SBarry Smith 
742dcf2fd19SBarry Smith   Input Parameters:
743420bcc1bSBarry Smith + ls           - `SNESLineSearch` object to monitor
744420bcc1bSBarry Smith . name         - the monitor type
745dcf2fd19SBarry Smith . help         - message indicating what monitoring is done
746dcf2fd19SBarry Smith . manual       - manual page for the monitor
747dcf2fd19SBarry Smith . monitor      - the monitor function
748f6dfbefdSBarry 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`
749dcf2fd19SBarry Smith 
750420bcc1bSBarry Smith   Calling sequence of `monitor`:
751420bcc1bSBarry Smith + ls - `SNESLineSearch` object being monitored
752420bcc1bSBarry Smith - vf - a `PetscViewerAndFormat` struct that provides the `PetscViewer` and `PetscViewerFormat` being used
753dcf2fd19SBarry Smith 
754420bcc1bSBarry Smith   Calling sequence of `monitorsetup`:
755420bcc1bSBarry Smith + ls - `SNESLineSearch` object being monitored
756420bcc1bSBarry Smith - vf - a `PetscViewerAndFormat` struct that provides the `PetscViewer` and `PetscViewerFormat` being used
757420bcc1bSBarry Smith 
758420bcc1bSBarry Smith   Level: advanced
759420bcc1bSBarry Smith 
760*648c30bcSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetMonitor()`, `PetscOptionsCreateViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
761db781477SPatrick Sanan           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
762e4094ef1SJacob Faibussowitsch           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
763db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
764c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
765db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
766db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
767dcf2fd19SBarry Smith @*/
768420bcc1bSBarry 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))
769d71ae5a4SJacob Faibussowitsch {
770dcf2fd19SBarry Smith   PetscViewer       viewer;
771dcf2fd19SBarry Smith   PetscViewerFormat format;
772dcf2fd19SBarry Smith   PetscBool         flg;
773dcf2fd19SBarry Smith 
774dcf2fd19SBarry Smith   PetscFunctionBegin;
775*648c30bcSBarry Smith   PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)ls), ((PetscObject)ls)->options, ((PetscObject)ls)->prefix, name, &viewer, &format, &flg));
776dcf2fd19SBarry Smith   if (flg) {
777d12e167eSBarry Smith     PetscViewerAndFormat *vf;
7789566063dSJacob Faibussowitsch     PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf));
779*648c30bcSBarry Smith     PetscCall(PetscViewerDestroy(&viewer));
7801baa6e33SBarry Smith     if (monitorsetup) PetscCall((*monitorsetup)(ls, vf));
7819566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchMonitorSet(ls, (PetscErrorCode(*)(SNESLineSearch, void *))monitor, vf, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy));
782dcf2fd19SBarry Smith   }
7833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
784dcf2fd19SBarry Smith }
785dcf2fd19SBarry Smith 
786f40b411bSPeter Brune /*@
787f1c6b773SPeter Brune   SNESLineSearchSetFromOptions - Sets options for the line search
788f40b411bSPeter Brune 
789c3339decSBarry Smith   Logically Collective
790f6dfbefdSBarry Smith 
791f899ff85SJose E. Roman   Input Parameter:
792420bcc1bSBarry Smith . linesearch - a `SNESLineSearch` line search context
793f40b411bSPeter Brune 
794f40b411bSPeter Brune   Options Database Keys:
7950b00b554SBarry Smith + -snes_linesearch_type <type>                                      - basic (or equivalently none), bt, l2, cp, nleqerr, shell
7963c7d6663SPeter Brune . -snes_linesearch_order <order>                                    - 1, 2, 3.  Most types only support certain orders (bt supports 2 or 3)
797f6dfbefdSBarry Smith . -snes_linesearch_norms                                            - Turn on/off the linesearch norms for the basic linesearch typem (`SNESLineSearchSetComputeNorms()`)
79871eef1aeSPeter Brune . -snes_linesearch_minlambda                                        - The minimum step length
7991a9b3a06SPeter Brune . -snes_linesearch_maxstep                                          - The maximum step size
8001a9b3a06SPeter Brune . -snes_linesearch_rtol                                             - Relative tolerance for iterative line searches
8011a9b3a06SPeter Brune . -snes_linesearch_atol                                             - Absolute tolerance for iterative line searches
8021a9b3a06SPeter Brune . -snes_linesearch_ltol                                             - Change in lambda tolerance for iterative line searches
8031a9b3a06SPeter Brune . -snes_linesearch_max_it                                           - The number of iterations for iterative line searches
804dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename]                              - Print progress of line searches
805dcf2fd19SBarry Smith . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine
8068e557f58SPeter Brune . -snes_linesearch_damping                                          - The linesearch damping parameter
807cd7522eaSPeter Brune . -snes_linesearch_keeplambda                                       - Keep the previous search length as the initial guess.
8081a9b3a06SPeter Brune . -snes_linesearch_precheck_picard                                  - Use precheck that speeds up convergence of picard method
809d8d34be6SBarry Smith - -snes_linesearch_precheck_picard_angle                            - Angle used in Picard precheck method
810f40b411bSPeter Brune 
811f40b411bSPeter Brune   Level: intermediate
812f40b411bSPeter Brune 
813420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchCreate()`, `SNESLineSearchSetOrder()`, `SNESLineSearchSetType()`, `SNESLineSearchSetTolerances()`, `SNESLineSearchSetDamping()`, `SNESLineSearchPreCheckPicard()`,
814db781477SPatrick Sanan           `SNESLineSearchType`, `SNESLineSearchSetComputeNorms()`
815f40b411bSPeter Brune @*/
816d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch)
817d71ae5a4SJacob Faibussowitsch {
8181a4f838cSPeter Brune   const char *deft = SNESLINESEARCHBASIC;
819bf7f4e0aSPeter Brune   char        type[256];
820bf7f4e0aSPeter Brune   PetscBool   flg, set;
821dcf2fd19SBarry Smith   PetscViewer viewer;
822bf388a1fSBarry Smith 
823bf7f4e0aSPeter Brune   PetscFunctionBegin;
8249566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchRegisterAll());
825bf7f4e0aSPeter Brune 
826d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)linesearch);
827f5af7f23SKarl Rupp   if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name;
8289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-snes_linesearch_type", "Linesearch type", "SNESLineSearchSetType", SNESLineSearchList, deft, type, 256, &flg));
829bf7f4e0aSPeter Brune   if (flg) {
8309566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetType(linesearch, type));
831bf7f4e0aSPeter Brune   } else if (!((PetscObject)linesearch)->type_name) {
8329566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetType(linesearch, deft));
833bf7f4e0aSPeter Brune   }
834bf7f4e0aSPeter Brune 
835*648c30bcSBarry Smith   PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)linesearch), ((PetscObject)linesearch)->options, ((PetscObject)linesearch)->prefix, "-snes_linesearch_monitor", &viewer, NULL, &set));
836*648c30bcSBarry Smith   if (set) PetscCall(SNESLineSearchSetDefaultMonitor(linesearch, viewer));
8379566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchMonitorSetFromOptions(linesearch, "-snes_linesearch_monitor_solution_update", "View correction at each iteration", "SNESLineSearchMonitorSolutionUpdate", SNESLineSearchMonitorSolutionUpdate, NULL));
838bf7f4e0aSPeter Brune 
8391a9b3a06SPeter Brune   /* tolerances */
8409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_minlambda", "Minimum step length", "SNESLineSearchSetTolerances", linesearch->steptol, &linesearch->steptol, NULL));
8419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_maxstep", "Maximum step size", "SNESLineSearchSetTolerances", linesearch->maxstep, &linesearch->maxstep, NULL));
8429566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_rtol", "Relative tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->rtol, &linesearch->rtol, NULL));
8439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_atol", "Absolute tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->atol, &linesearch->atol, NULL));
8449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_ltol", "Change in lambda tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->ltol, &linesearch->ltol, NULL));
8459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_linesearch_max_it", "Maximum iterations for iterative line searches", "SNESLineSearchSetTolerances", linesearch->max_its, &linesearch->max_its, NULL));
8467a35526eSPeter Brune 
8471a9b3a06SPeter Brune   /* damping parameters */
8489566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_damping", "Line search damping and initial step guess", "SNESLineSearchSetDamping", linesearch->damping, &linesearch->damping, NULL));
8491a9b3a06SPeter Brune 
8509566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_linesearch_keeplambda", "Use previous lambda as damping", "SNESLineSearchSetKeepLambda", linesearch->keeplambda, &linesearch->keeplambda, NULL));
8511a9b3a06SPeter Brune 
8521a9b3a06SPeter Brune   /* precheck */
8539566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_linesearch_precheck_picard", "Use a correction that sometimes improves convergence of Picard iteration", "SNESLineSearchPreCheckPicard", flg, &flg, &set));
8547a35526eSPeter Brune   if (set) {
8557a35526eSPeter Brune     if (flg) {
8567a35526eSPeter Brune       linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */
857f5af7f23SKarl Rupp 
858d0609cedSBarry 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));
8599566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchSetPreCheck(linesearch, SNESLineSearchPreCheckPicard, &linesearch->precheck_picard_angle));
8607a35526eSPeter Brune     } else {
8619566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchSetPreCheck(linesearch, NULL, NULL));
8627a35526eSPeter Brune     }
8637a35526eSPeter Brune   }
8649566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_linesearch_order", "Order of approximation used in the line search", "SNESLineSearchSetOrder", linesearch->order, &linesearch->order, NULL));
8659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_linesearch_norms", "Compute final norms in line search", "SNESLineSearchSetComputeNorms", linesearch->norms, &linesearch->norms, NULL));
8667a35526eSPeter Brune 
867dbbe0bcdSBarry Smith   PetscTryTypeMethod(linesearch, setfromoptions, PetscOptionsObject);
8685a9b6599SPeter Brune 
869dbbe0bcdSBarry Smith   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)linesearch, PetscOptionsObject));
870d0609cedSBarry Smith   PetscOptionsEnd();
8713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
872bf7f4e0aSPeter Brune }
873bf7f4e0aSPeter Brune 
874f40b411bSPeter Brune /*@
875f190f2fcSBarry Smith   SNESLineSearchView - Prints useful information about the line search
876f40b411bSPeter Brune 
87720f4b53cSBarry Smith   Logically Collective
87820f4b53cSBarry Smith 
879f40b411bSPeter Brune   Input Parameters:
8802fe279fdSBarry Smith + linesearch - line search context
881420bcc1bSBarry Smith - viewer     - the `PetscViewer` to display the line search information to
882f40b411bSPeter Brune 
883f40b411bSPeter Brune   Level: intermediate
884f40b411bSPeter Brune 
885420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `PetscViewer`, `SNESLineSearchCreate()`
886f40b411bSPeter Brune @*/
887d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer)
888d71ae5a4SJacob Faibussowitsch {
8897f1410a3SPeter Brune   PetscBool iascii;
890bf388a1fSBarry Smith 
891bf7f4e0aSPeter Brune   PetscFunctionBegin;
8927f1410a3SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
89348a46eb9SPierre Jolivet   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch), &viewer));
8947f1410a3SPeter Brune   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
8957f1410a3SPeter Brune   PetscCheckSameComm(linesearch, 1, viewer, 2);
896f40b411bSPeter Brune 
8979566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
8987f1410a3SPeter Brune   if (iascii) {
8999566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)linesearch, viewer));
9009566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
901dbbe0bcdSBarry Smith     PetscTryTypeMethod(linesearch, view, viewer);
9029566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
9039566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep, (double)linesearch->steptol));
9049566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol, (double)linesearch->atol, (double)linesearch->ltol));
90563a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  maximum iterations=%" PetscInt_FMT "\n", linesearch->max_its));
9066b2b7091SBarry Smith     if (linesearch->ops->precheck) {
9076b2b7091SBarry Smith       if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) {
90863a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  using precheck step to speed up Picard convergence\n"));
9097f1410a3SPeter Brune       } else {
91063a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  using user-defined precheck step\n"));
9117f1410a3SPeter Brune       }
9127f1410a3SPeter Brune     }
91348a46eb9SPierre Jolivet     if (linesearch->ops->postcheck) PetscCall(PetscViewerASCIIPrintf(viewer, "  using user-defined postcheck step\n"));
9147f1410a3SPeter Brune   }
9153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
916bf7f4e0aSPeter Brune }
917bf7f4e0aSPeter Brune 
918cc4c1da9SBarry Smith /*@
919420bcc1bSBarry Smith   SNESLineSearchGetType - Gets the `SNESLinesearchType` of a `SNESLineSearch`
920a80ff896SJed Brown 
921c3339decSBarry Smith   Logically Collective
922a80ff896SJed Brown 
9232fe279fdSBarry Smith   Input Parameter:
924420bcc1bSBarry Smith . linesearch - the line search context
925a80ff896SJed Brown 
9262fe279fdSBarry Smith   Output Parameter:
9272fe279fdSBarry Smith . type - The type of line search, or `NULL` if not set
928a80ff896SJed Brown 
929a80ff896SJed Brown   Level: intermediate
930a80ff896SJed Brown 
931420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchSetType()`
932a80ff896SJed Brown @*/
933d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetType(SNESLineSearch linesearch, SNESLineSearchType *type)
934d71ae5a4SJacob Faibussowitsch {
935a80ff896SJed Brown   PetscFunctionBegin;
936a80ff896SJed Brown   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
9374f572ea9SToby Isaac   PetscAssertPointer(type, 2);
938a80ff896SJed Brown   *type = ((PetscObject)linesearch)->type_name;
9393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
940a80ff896SJed Brown }
941a80ff896SJed Brown 
942cc4c1da9SBarry Smith /*@
943420bcc1bSBarry Smith   SNESLineSearchSetType - Sets the `SNESLinesearchType` of a `SNESLineSearch`
944f40b411bSPeter Brune 
945c3339decSBarry Smith   Logically Collective
946f190f2fcSBarry Smith 
947f40b411bSPeter Brune   Input Parameters:
948420bcc1bSBarry Smith + linesearch - the line search context
949ceaaa498SBarry Smith - type       - The type of line search to be used, see `SNESLineSearchType`
9501fe24845SBarry Smith 
9513c7db156SBarry Smith   Options Database Key:
9520b00b554SBarry Smith . -snes_linesearch_type <type> - basic (or equivalently none), bt, l2, cp, nleqerr, shell
953cd7522eaSPeter Brune 
954f40b411bSPeter Brune   Level: intermediate
955f40b411bSPeter Brune 
956420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchGetType()`
957f40b411bSPeter Brune @*/
958d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type)
959d71ae5a4SJacob Faibussowitsch {
960bf7f4e0aSPeter Brune   PetscBool match;
9615f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(SNESLineSearch);
962bf7f4e0aSPeter Brune 
963bf7f4e0aSPeter Brune   PetscFunctionBegin;
964f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
9654f572ea9SToby Isaac   PetscAssertPointer(type, 2);
966bf7f4e0aSPeter Brune 
9679566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)linesearch, type, &match));
9683ba16761SJacob Faibussowitsch   if (match) PetscFunctionReturn(PETSC_SUCCESS);
969bf7f4e0aSPeter Brune 
9709566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(SNESLineSearchList, type, &r));
9716adde796SStefano Zampini   PetscCheck(r, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested Line Search type %s", type);
972bf7f4e0aSPeter Brune   /* Destroy the previous private line search context */
973213acdd3SPierre Jolivet   PetscTryTypeMethod(linesearch, destroy);
9740298fd71SBarry Smith   linesearch->ops->destroy = NULL;
975f1c6b773SPeter Brune   /* Reinitialize function pointers in SNESLineSearchOps structure */
9769e5d0892SLisandro Dalcin   linesearch->ops->apply          = NULL;
9779e5d0892SLisandro Dalcin   linesearch->ops->view           = NULL;
9789e5d0892SLisandro Dalcin   linesearch->ops->setfromoptions = NULL;
9799e5d0892SLisandro Dalcin   linesearch->ops->destroy        = NULL;
980bf7f4e0aSPeter Brune 
9819566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)linesearch, type));
9829566063dSJacob Faibussowitsch   PetscCall((*r)(linesearch));
9833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
984bf7f4e0aSPeter Brune }
985bf7f4e0aSPeter Brune 
986f40b411bSPeter Brune /*@
987f6dfbefdSBarry Smith   SNESLineSearchSetSNES - Sets the `SNES` for the linesearch for function evaluation.
988f40b411bSPeter Brune 
989f40b411bSPeter Brune   Input Parameters:
990420bcc1bSBarry Smith + linesearch - the line search context
991420bcc1bSBarry Smith - snes       - The `SNES` instance
992f40b411bSPeter Brune 
99378bcb3b5SPeter Brune   Level: developer
99478bcb3b5SPeter Brune 
995f6dfbefdSBarry Smith   Note:
996f190f2fcSBarry Smith   This happens automatically when the line search is obtained/created with
997f6dfbefdSBarry Smith   `SNESGetLineSearch()`.  This routine is therefore mainly called within `SNES`
99878bcb3b5SPeter Brune   implementations.
999f40b411bSPeter Brune 
1000420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetSNES()`, `SNESLineSearchSetVecs()`
1001f40b411bSPeter Brune @*/
1002d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes)
1003d71ae5a4SJacob Faibussowitsch {
1004bf7f4e0aSPeter Brune   PetscFunctionBegin;
1005f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1006bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 2);
1007bf7f4e0aSPeter Brune   linesearch->snes = snes;
10083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1009bf7f4e0aSPeter Brune }
1010bf7f4e0aSPeter Brune 
1011f40b411bSPeter Brune /*@
1012f6dfbefdSBarry Smith   SNESLineSearchGetSNES - Gets the `SNES` instance associated with the line search.
1013f6dfbefdSBarry Smith 
1014f6dfbefdSBarry Smith   Not Collective
1015f40b411bSPeter Brune 
10162fe279fdSBarry Smith   Input Parameter:
1017420bcc1bSBarry Smith . linesearch - the line search context
1018f40b411bSPeter Brune 
10192fe279fdSBarry Smith   Output Parameter:
10202fe279fdSBarry Smith . snes - The `SNES` instance
1021f40b411bSPeter Brune 
10228141a3b9SPeter Brune   Level: developer
1023f40b411bSPeter Brune 
1024420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESType`, `SNESLineSearchSetVecs()`
1025f40b411bSPeter Brune @*/
1026d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes)
1027d71ae5a4SJacob Faibussowitsch {
1028bf7f4e0aSPeter Brune   PetscFunctionBegin;
1029f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
10304f572ea9SToby Isaac   PetscAssertPointer(snes, 2);
1031bf7f4e0aSPeter Brune   *snes = linesearch->snes;
10323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1033bf7f4e0aSPeter Brune }
1034bf7f4e0aSPeter Brune 
1035f40b411bSPeter Brune /*@
1036420bcc1bSBarry Smith   SNESLineSearchGetLambda - Gets the last line search steplength used
1037f40b411bSPeter Brune 
1038f6dfbefdSBarry Smith   Not Collective
1039f6dfbefdSBarry Smith 
10402fe279fdSBarry Smith   Input Parameter:
1041420bcc1bSBarry Smith . linesearch - the line search context
1042f40b411bSPeter Brune 
10432fe279fdSBarry Smith   Output Parameter:
1044f6dfbefdSBarry Smith . lambda - The last steplength computed during `SNESLineSearchApply()`
1045f40b411bSPeter Brune 
104678bcb3b5SPeter Brune   Level: advanced
104778bcb3b5SPeter Brune 
1048f6dfbefdSBarry Smith   Note:
10498e557f58SPeter Brune   This is useful in methods where the solver is ill-scaled and
105078bcb3b5SPeter Brune   requires some adaptive notion of the difference in scale between the
1051f6dfbefdSBarry Smith   solution and the function.  For instance, `SNESQN` may be scaled by the
105278bcb3b5SPeter Brune   line search lambda using the argument -snes_qn_scaling ls.
105378bcb3b5SPeter Brune 
1054420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetLambda()`, `SNESLineSearchGetDamping()`, `SNESLineSearchApply()`
1055f40b411bSPeter Brune @*/
1056d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch linesearch, PetscReal *lambda)
1057d71ae5a4SJacob Faibussowitsch {
10586a388c36SPeter Brune   PetscFunctionBegin;
1059f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
10604f572ea9SToby Isaac   PetscAssertPointer(lambda, 2);
10616a388c36SPeter Brune   *lambda = linesearch->lambda;
10623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10636a388c36SPeter Brune }
10646a388c36SPeter Brune 
1065f40b411bSPeter Brune /*@
1066f6dfbefdSBarry Smith   SNESLineSearchSetLambda - Sets the line search steplength
1067f40b411bSPeter Brune 
1068f40b411bSPeter Brune   Input Parameters:
10698e557f58SPeter Brune + linesearch - line search context
1070420bcc1bSBarry Smith - lambda     - The steplength to use
1071f40b411bSPeter Brune 
107220f4b53cSBarry Smith   Level: advanced
107320f4b53cSBarry Smith 
1074f6dfbefdSBarry Smith   Note:
1075f6dfbefdSBarry Smith   This routine is typically used within implementations of `SNESLineSearchApply()`
1076f6dfbefdSBarry Smith   to set the final steplength.  This routine (and `SNESLineSearchGetLambda()`) were
1077cd7522eaSPeter Brune   added in order to facilitate Quasi-Newton methods that use the previous steplength
1078cd7522eaSPeter Brune   as an inner scaling parameter.
1079cd7522eaSPeter Brune 
1080420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetLambda()`
1081f40b411bSPeter Brune @*/
1082d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
1083d71ae5a4SJacob Faibussowitsch {
10846a388c36SPeter Brune   PetscFunctionBegin;
1085f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
10866a388c36SPeter Brune   linesearch->lambda = lambda;
10873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10886a388c36SPeter Brune }
10896a388c36SPeter Brune 
1090f40b411bSPeter Brune /*@
1091ceaaa498SBarry Smith   SNESLineSearchGetTolerances - Gets the tolerances for the line search.
1092f40b411bSPeter Brune 
1093f6dfbefdSBarry Smith   Not Collective
1094f6dfbefdSBarry Smith 
1095f899ff85SJose E. Roman   Input Parameter:
1096420bcc1bSBarry Smith . linesearch - the line search context
1097f40b411bSPeter Brune 
1098f40b411bSPeter Brune   Output Parameters:
1099b13b64c2SBarry Smith + steptol - The minimum steplength
11006cc8e53bSPeter Brune . maxstep - The maximum steplength
1101516fe3c3SPeter Brune . rtol    - The relative tolerance for iterative line searches
1102516fe3c3SPeter Brune . atol    - The absolute tolerance for iterative line searches
1103516fe3c3SPeter Brune . ltol    - The change in lambda tolerance for iterative line searches
1104e4094ef1SJacob Faibussowitsch - max_its - The maximum number of iterations of the line search
1105f40b411bSPeter Brune 
110678bcb3b5SPeter Brune   Level: intermediate
110778bcb3b5SPeter Brune 
1108f6dfbefdSBarry Smith   Note:
110978bcb3b5SPeter Brune   Different line searches may implement these parameters slightly differently as
11103c7d6663SPeter Brune   the type requires.
1111516fe3c3SPeter Brune 
1112420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetTolerances()`
1113f40b411bSPeter Brune @*/
1114d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch linesearch, PetscReal *steptol, PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
1115d71ae5a4SJacob Faibussowitsch {
11166a388c36SPeter Brune   PetscFunctionBegin;
1117f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1118516fe3c3SPeter Brune   if (steptol) {
11194f572ea9SToby Isaac     PetscAssertPointer(steptol, 2);
11206a388c36SPeter Brune     *steptol = linesearch->steptol;
1121516fe3c3SPeter Brune   }
1122516fe3c3SPeter Brune   if (maxstep) {
11234f572ea9SToby Isaac     PetscAssertPointer(maxstep, 3);
1124b13b64c2SBarry Smith     *maxstep = linesearch->maxstep;
1125516fe3c3SPeter Brune   }
1126516fe3c3SPeter Brune   if (rtol) {
11274f572ea9SToby Isaac     PetscAssertPointer(rtol, 4);
1128516fe3c3SPeter Brune     *rtol = linesearch->rtol;
1129516fe3c3SPeter Brune   }
1130516fe3c3SPeter Brune   if (atol) {
11314f572ea9SToby Isaac     PetscAssertPointer(atol, 5);
1132516fe3c3SPeter Brune     *atol = linesearch->atol;
1133516fe3c3SPeter Brune   }
1134516fe3c3SPeter Brune   if (ltol) {
11354f572ea9SToby Isaac     PetscAssertPointer(ltol, 6);
1136516fe3c3SPeter Brune     *ltol = linesearch->ltol;
1137516fe3c3SPeter Brune   }
1138516fe3c3SPeter Brune   if (max_its) {
11394f572ea9SToby Isaac     PetscAssertPointer(max_its, 7);
1140516fe3c3SPeter Brune     *max_its = linesearch->max_its;
1141516fe3c3SPeter Brune   }
11423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11436a388c36SPeter Brune }
11446a388c36SPeter Brune 
1145f40b411bSPeter Brune /*@
1146ceaaa498SBarry Smith   SNESLineSearchSetTolerances -  Sets the tolerances for the linesearch.
1147f40b411bSPeter Brune 
1148f6dfbefdSBarry Smith   Collective
1149f6dfbefdSBarry Smith 
1150f40b411bSPeter Brune   Input Parameters:
1151420bcc1bSBarry Smith + linesearch - the line search context
1152516fe3c3SPeter Brune . steptol    - The minimum steplength
11536cc8e53bSPeter Brune . maxstep    - The maximum steplength
1154516fe3c3SPeter Brune . rtol       - The relative tolerance for iterative line searches
1155516fe3c3SPeter Brune . atol       - The absolute tolerance for iterative line searches
1156516fe3c3SPeter Brune . ltol       - The change in lambda tolerance for iterative line searches
1157420bcc1bSBarry Smith - max_it     - The maximum number of iterations of the line search
1158420bcc1bSBarry Smith 
1159420bcc1bSBarry Smith   Options Database Keys:
1160420bcc1bSBarry Smith + -snes_linesearch_minlambda - The minimum step length
1161420bcc1bSBarry Smith . -snes_linesearch_maxstep   - The maximum step size
1162420bcc1bSBarry Smith . -snes_linesearch_rtol      - Relative tolerance for iterative line searches
1163420bcc1bSBarry Smith . -snes_linesearch_atol      - Absolute tolerance for iterative line searches
1164420bcc1bSBarry Smith . -snes_linesearch_ltol      - Change in lambda tolerance for iterative line searches
1165420bcc1bSBarry Smith - -snes_linesearch_max_it    - The number of iterations for iterative line searches
1166f40b411bSPeter Brune 
116720f4b53cSBarry Smith   Level: intermediate
116820f4b53cSBarry Smith 
1169f6dfbefdSBarry Smith   Note:
1170420bcc1bSBarry Smith   The user may choose to not set any of the tolerances using `PETSC_DEFAULT` in place of an argument.
1171f40b411bSPeter Brune 
1172420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetTolerances()`
1173f40b411bSPeter Brune @*/
1174420bcc1bSBarry Smith PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch linesearch, PetscReal steptol, PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_it)
1175d71ae5a4SJacob Faibussowitsch {
11766a388c36SPeter Brune   PetscFunctionBegin;
1177f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1178d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, steptol, 2);
1179d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, maxstep, 3);
1180d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, rtol, 4);
1181d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, atol, 5);
1182d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, ltol, 6);
1183420bcc1bSBarry Smith   PetscValidLogicalCollectiveInt(linesearch, max_it, 7);
1184d3952184SSatish Balay 
118513bcc0bdSJacob Faibussowitsch   if (steptol != (PetscReal)PETSC_DEFAULT) {
11865f80ce2aSJacob Faibussowitsch     PetscCheck(steptol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Minimum step length %14.12e must be non-negative", (double)steptol);
11876a388c36SPeter Brune     linesearch->steptol = steptol;
1188d3952184SSatish Balay   }
1189d3952184SSatish Balay 
119013bcc0bdSJacob Faibussowitsch   if (maxstep != (PetscReal)PETSC_DEFAULT) {
11915f80ce2aSJacob Faibussowitsch     PetscCheck(maxstep >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum step length %14.12e must be non-negative", (double)maxstep);
1192516fe3c3SPeter Brune     linesearch->maxstep = maxstep;
1193d3952184SSatish Balay   }
1194d3952184SSatish Balay 
119513bcc0bdSJacob Faibussowitsch   if (rtol != (PetscReal)PETSC_DEFAULT) {
11962061ca32SJunchao 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);
1197516fe3c3SPeter Brune     linesearch->rtol = rtol;
1198d3952184SSatish Balay   }
1199d3952184SSatish Balay 
120013bcc0bdSJacob Faibussowitsch   if (atol != (PetscReal)PETSC_DEFAULT) {
12015f80ce2aSJacob Faibussowitsch     PetscCheck(atol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Absolute tolerance %14.12e must be non-negative", (double)atol);
1202516fe3c3SPeter Brune     linesearch->atol = atol;
1203d3952184SSatish Balay   }
1204d3952184SSatish Balay 
120513bcc0bdSJacob Faibussowitsch   if (ltol != (PetscReal)PETSC_DEFAULT) {
12065f80ce2aSJacob Faibussowitsch     PetscCheck(ltol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Lambda tolerance %14.12e must be non-negative", (double)ltol);
1207516fe3c3SPeter Brune     linesearch->ltol = ltol;
1208d3952184SSatish Balay   }
1209d3952184SSatish Balay 
1210420bcc1bSBarry Smith   if (max_it != PETSC_DEFAULT) {
1211420bcc1bSBarry Smith     PetscCheck(max_it >= 0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of iterations %" PetscInt_FMT " must be non-negative", max_it);
1212420bcc1bSBarry Smith     linesearch->max_its = max_it;
1213d3952184SSatish Balay   }
12143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12156a388c36SPeter Brune }
12166a388c36SPeter Brune 
1217f40b411bSPeter Brune /*@
1218f1c6b773SPeter Brune   SNESLineSearchGetDamping - Gets the line search damping parameter.
1219f40b411bSPeter Brune 
12202fe279fdSBarry Smith   Input Parameter:
1221420bcc1bSBarry Smith . linesearch - the line search context
1222f40b411bSPeter Brune 
12232fe279fdSBarry Smith   Output Parameter:
12248e557f58SPeter Brune . damping - The damping parameter
1225f40b411bSPeter Brune 
122678bcb3b5SPeter Brune   Level: advanced
1227f40b411bSPeter Brune 
1228420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearchGetStepTolerance()`, `SNESQN`
1229f40b411bSPeter Brune @*/
1230d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch linesearch, PetscReal *damping)
1231d71ae5a4SJacob Faibussowitsch {
12326a388c36SPeter Brune   PetscFunctionBegin;
1233f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
12344f572ea9SToby Isaac   PetscAssertPointer(damping, 2);
12356a388c36SPeter Brune   *damping = linesearch->damping;
12363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12376a388c36SPeter Brune }
12386a388c36SPeter Brune 
1239f40b411bSPeter Brune /*@
1240fd292e60Sprj-   SNESLineSearchSetDamping - Sets the line search damping parameter.
1241f40b411bSPeter Brune 
1242f40b411bSPeter Brune   Input Parameters:
1243420bcc1bSBarry Smith + linesearch - the line search context
124403fd4120SBarry Smith - damping    - The damping parameter
1245f40b411bSPeter Brune 
1246f6dfbefdSBarry Smith   Options Database Key:
1247f6dfbefdSBarry Smith . -snes_linesearch_damping <damping> - the damping value
1248f6dfbefdSBarry Smith 
1249f40b411bSPeter Brune   Level: intermediate
1250f40b411bSPeter Brune 
1251f6dfbefdSBarry Smith   Note:
1252f6dfbefdSBarry Smith   The `SNESLINESEARCHNONE` line search merely takes the update step scaled by the damping parameter.
1253420bcc1bSBarry Smith   The use of the damping parameter in the `SNESLINESEARCHL2` and `SNESLINESEARCHCP` line searches is much more subtle;
125478bcb3b5SPeter Brune   it is used as a starting point in calculating the secant step. However, the eventual
1255420bcc1bSBarry Smith   step may be of greater length than the damping parameter.  In the `SNESLINESEARCHBT` line search it is
1256420bcc1bSBarry Smith   used as the maximum possible step length, as the `SNESLINESEARCHBT` line search only backtracks.
1257cd7522eaSPeter Brune 
1258420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetDamping()`
1259f40b411bSPeter Brune @*/
1260d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch linesearch, PetscReal damping)
1261d71ae5a4SJacob Faibussowitsch {
12626a388c36SPeter Brune   PetscFunctionBegin;
1263f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
12646a388c36SPeter Brune   linesearch->damping = damping;
12653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12666a388c36SPeter Brune }
12676a388c36SPeter Brune 
126859405d5eSPeter Brune /*@
126959405d5eSPeter Brune   SNESLineSearchGetOrder - Gets the line search approximation order.
127059405d5eSPeter Brune 
1271f6dfbefdSBarry Smith   Input Parameter:
1272420bcc1bSBarry Smith . linesearch - the line search context
127359405d5eSPeter Brune 
1274f6dfbefdSBarry Smith   Output Parameter:
12758e557f58SPeter Brune . order - The order
127659405d5eSPeter Brune 
127759405d5eSPeter Brune   Level: intermediate
127859405d5eSPeter Brune 
1279420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetOrder()`
128059405d5eSPeter Brune @*/
1281d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch linesearch, PetscInt *order)
1282d71ae5a4SJacob Faibussowitsch {
128359405d5eSPeter Brune   PetscFunctionBegin;
128459405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
12854f572ea9SToby Isaac   PetscAssertPointer(order, 2);
128659405d5eSPeter Brune   *order = linesearch->order;
12873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
128859405d5eSPeter Brune }
128959405d5eSPeter Brune 
129059405d5eSPeter Brune /*@
12911f8196a2SJed Brown   SNESLineSearchSetOrder - Sets the maximum order of the polynomial fit used in the line search
129259405d5eSPeter Brune 
129359405d5eSPeter Brune   Input Parameters:
1294420bcc1bSBarry Smith + linesearch - the line search context
1295ceaaa498SBarry Smith - order      - The order
129659405d5eSPeter Brune 
129759405d5eSPeter Brune   Level: intermediate
129859405d5eSPeter Brune 
1299ceaaa498SBarry Smith   Values for `order`\:
1300f6dfbefdSBarry Smith +  1 or `SNES_LINESEARCH_ORDER_LINEAR` - linear order
1301f6dfbefdSBarry Smith .  2 or `SNES_LINESEARCH_ORDER_QUADRATIC` - quadratic order
1302f6dfbefdSBarry Smith -  3 or `SNES_LINESEARCH_ORDER_CUBIC` - cubic order
130378bcb3b5SPeter Brune 
1304420bcc1bSBarry Smith   Options Database Key:
1305420bcc1bSBarry Smith . -snes_linesearch_order <order> - 1, 2, 3.  Most types only support certain orders (`SNESLINESEARCHBT` supports 2 or 3)
1306420bcc1bSBarry Smith 
1307ceaaa498SBarry Smith   Note:
1308ceaaa498SBarry Smith   These orders are supported by `SNESLINESEARCHBT` and `SNESLINESEARCHCP`
130959405d5eSPeter Brune 
1310420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetOrder()`, `SNESLineSearchSetDamping()`
131159405d5eSPeter Brune @*/
1312d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch linesearch, PetscInt order)
1313d71ae5a4SJacob Faibussowitsch {
131459405d5eSPeter Brune   PetscFunctionBegin;
131559405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
131659405d5eSPeter Brune   linesearch->order = order;
13173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
131859405d5eSPeter Brune }
131959405d5eSPeter Brune 
1320f40b411bSPeter Brune /*@
1321420bcc1bSBarry Smith   SNESLineSearchGetNorms - Gets the norms for the current solution `X`, the current update `Y`, and the current function value `F`.
1322f40b411bSPeter Brune 
1323f6dfbefdSBarry Smith   Not Collective
1324f6dfbefdSBarry Smith 
1325f899ff85SJose E. Roman   Input Parameter:
1326420bcc1bSBarry Smith . linesearch - the line search context
1327f40b411bSPeter Brune 
1328f40b411bSPeter Brune   Output Parameters:
1329f40b411bSPeter Brune + xnorm - The norm of the current solution
1330a68bbae5SBarry Smith . fnorm - The norm of the current function, this is the `norm(function(X))` where `X` is the current solution.
13316b095a85SStefano Zampini - ynorm - The norm of the current update (after scaling by the linesearch computed lambda)
1332f40b411bSPeter Brune 
133378bcb3b5SPeter Brune   Level: developer
1334f40b411bSPeter Brune 
1335420bcc1bSBarry Smith   Notes:
1336420bcc1bSBarry Smith   Some values may not be up-to-date at particular points in the code.
1337a68bbae5SBarry Smith 
1338a68bbae5SBarry Smith   This, in combination with `SNESLineSearchSetNorms()`, allow the line search and the `SNESSolve_XXX()` to share
1339a68bbae5SBarry Smith   computed values.
1340a68bbae5SBarry Smith 
1341420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetNorms()` `SNESLineSearchGetVecs()`
1342f40b411bSPeter Brune @*/
1343d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal *xnorm, PetscReal *fnorm, PetscReal *ynorm)
1344d71ae5a4SJacob Faibussowitsch {
1345bf7f4e0aSPeter Brune   PetscFunctionBegin;
1346f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1347f5af7f23SKarl Rupp   if (xnorm) *xnorm = linesearch->xnorm;
1348f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
1349f5af7f23SKarl Rupp   if (ynorm) *ynorm = linesearch->ynorm;
13503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1351bf7f4e0aSPeter Brune }
1352bf7f4e0aSPeter Brune 
1353f40b411bSPeter Brune /*@
13541bd63e3eSSatish Balay   SNESLineSearchSetNorms - Sets the computed norms for the current solution `X`, the current update `Y`, and the current function value `F`.
1355f40b411bSPeter Brune 
1356c3339decSBarry Smith   Collective
1357f6dfbefdSBarry Smith 
1358f40b411bSPeter Brune   Input Parameters:
1359420bcc1bSBarry Smith + linesearch - the line search context
1360f40b411bSPeter Brune . xnorm      - The norm of the current solution
1361a68bbae5SBarry Smith . fnorm      - The norm of the current function, this is the `norm(function(X))` where `X` is the current solution
13626b095a85SStefano Zampini - ynorm      - The norm of the current update (after scaling by the linesearch computed lambda)
1363f40b411bSPeter Brune 
1364f6dfbefdSBarry Smith   Level: developer
1365f40b411bSPeter Brune 
1366420bcc1bSBarry Smith   Note:
1367420bcc1bSBarry Smith   This is called by the line search routines to store the values they have just computed
1368420bcc1bSBarry Smith 
1369420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()`
1370f40b411bSPeter Brune @*/
1371d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
1372d71ae5a4SJacob Faibussowitsch {
13736a388c36SPeter Brune   PetscFunctionBegin;
1374f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
13756a388c36SPeter Brune   linesearch->xnorm = xnorm;
13766a388c36SPeter Brune   linesearch->fnorm = fnorm;
13776a388c36SPeter Brune   linesearch->ynorm = ynorm;
13783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13796a388c36SPeter Brune }
13806a388c36SPeter Brune 
1381f40b411bSPeter Brune /*@
1382420bcc1bSBarry Smith   SNESLineSearchComputeNorms - Explicitly computes the norms of the current solution `X`, the current update `Y`, and the current function value `F`.
1383f40b411bSPeter Brune 
13842fe279fdSBarry Smith   Input Parameter:
1385420bcc1bSBarry Smith . linesearch - the line search context
1386f40b411bSPeter Brune 
138720f4b53cSBarry Smith   Options Database Key:
13888e557f58SPeter Brune . -snes_linesearch_norms - turn norm computation on or off
1389f40b411bSPeter Brune 
1390f40b411bSPeter Brune   Level: intermediate
1391f40b411bSPeter Brune 
1392420bcc1bSBarry Smith   Developer Note:
1393420bcc1bSBarry Smith   The options database key is misnamed. It should be -snes_linesearch_compute_norms
1394420bcc1bSBarry Smith 
1395420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms`, `SNESLineSearchSetNorms()`, `SNESLineSearchSetComputeNorms()`
1396f40b411bSPeter Brune @*/
1397d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
1398d71ae5a4SJacob Faibussowitsch {
13999bd66eb0SPeter Brune   SNES snes;
1400bf388a1fSBarry Smith 
14016a388c36SPeter Brune   PetscFunctionBegin;
14026a388c36SPeter Brune   if (linesearch->norms) {
14039bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
14049566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
14059566063dSJacob Faibussowitsch       PetscCall(VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
14069566063dSJacob Faibussowitsch       PetscCall(VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm));
14079566063dSJacob Faibussowitsch       PetscCall((*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm));
14089bd66eb0SPeter Brune     } else {
14099566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(linesearch->vec_func, NORM_2, &linesearch->fnorm));
14109566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
14119566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm));
14129566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(linesearch->vec_func, NORM_2, &linesearch->fnorm));
14139566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
14149566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(linesearch->vec_update, NORM_2, &linesearch->ynorm));
14156a388c36SPeter Brune     }
14169bd66eb0SPeter Brune   }
14173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14186a388c36SPeter Brune }
14196a388c36SPeter Brune 
14206f263ca3SPeter Brune /*@
14216f263ca3SPeter Brune   SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search.
14226f263ca3SPeter Brune 
14236f263ca3SPeter Brune   Input Parameters:
1424420bcc1bSBarry Smith + linesearch - the line search context
142578bcb3b5SPeter Brune - flg        - indicates whether or not to compute norms
14266f263ca3SPeter Brune 
142720f4b53cSBarry Smith   Options Database Key:
1428420bcc1bSBarry Smith . -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic (none) `SNESLINESEARCHBASIC` line search
14296f263ca3SPeter Brune 
143020f4b53cSBarry Smith   Level: intermediate
143120f4b53cSBarry Smith 
1432f6dfbefdSBarry Smith   Note:
1433f6dfbefdSBarry 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.
14346f263ca3SPeter Brune 
1435420bcc1bSBarry Smith   Developer Note:
1436420bcc1bSBarry Smith   The options database key is misnamed. It should be -snes_linesearch_compute_norms
1437420bcc1bSBarry Smith 
1438420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetNorms()`, `SNESLineSearchComputeNorms()`, `SNESLINESEARCHBASIC`
14396f263ca3SPeter Brune @*/
1440d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg)
1441d71ae5a4SJacob Faibussowitsch {
14426f263ca3SPeter Brune   PetscFunctionBegin;
14436f263ca3SPeter Brune   linesearch->norms = flg;
14443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14456f263ca3SPeter Brune }
14466f263ca3SPeter Brune 
1447f40b411bSPeter Brune /*@
1448f6dfbefdSBarry Smith   SNESLineSearchGetVecs - Gets the vectors from the `SNESLineSearch` context
1449f6dfbefdSBarry Smith 
14508f14a041SBarry Smith   Not Collective but the vectors are parallel
1451f40b411bSPeter Brune 
1452f899ff85SJose E. Roman   Input Parameter:
1453420bcc1bSBarry Smith . linesearch - the line search context
1454f40b411bSPeter Brune 
1455f40b411bSPeter Brune   Output Parameters:
14566232e825SPeter Brune + X - Solution vector
14576232e825SPeter Brune . F - Function vector
14586232e825SPeter Brune . Y - Search direction vector
14596232e825SPeter Brune . W - Solution work vector
14606232e825SPeter Brune - G - Function work vector
14616232e825SPeter Brune 
146220f4b53cSBarry Smith   Level: advanced
146320f4b53cSBarry Smith 
14647bba9028SPeter Brune   Notes:
146520f4b53cSBarry Smith   At the beginning of a line search application, `X` should contain a
146620f4b53cSBarry Smith   solution and the vector `F` the function computed at `X`.  At the end of the
146720f4b53cSBarry Smith   line search application, `X` should contain the new solution, and `F` the
14686232e825SPeter Brune   function evaluated at the new solution.
1469f40b411bSPeter Brune 
1470f6dfbefdSBarry Smith   These vectors are owned by the `SNESLineSearch` and should not be destroyed by the caller
14712a7a6963SBarry Smith 
1472420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()`
1473f40b411bSPeter Brune @*/
1474d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch, Vec *X, Vec *F, Vec *Y, Vec *W, Vec *G)
1475d71ae5a4SJacob Faibussowitsch {
14766a388c36SPeter Brune   PetscFunctionBegin;
1477f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
14786a388c36SPeter Brune   if (X) {
14794f572ea9SToby Isaac     PetscAssertPointer(X, 2);
14806a388c36SPeter Brune     *X = linesearch->vec_sol;
14816a388c36SPeter Brune   }
14826a388c36SPeter Brune   if (F) {
14834f572ea9SToby Isaac     PetscAssertPointer(F, 3);
14846a388c36SPeter Brune     *F = linesearch->vec_func;
14856a388c36SPeter Brune   }
14866a388c36SPeter Brune   if (Y) {
14874f572ea9SToby Isaac     PetscAssertPointer(Y, 4);
14886a388c36SPeter Brune     *Y = linesearch->vec_update;
14896a388c36SPeter Brune   }
14906a388c36SPeter Brune   if (W) {
14914f572ea9SToby Isaac     PetscAssertPointer(W, 5);
14926a388c36SPeter Brune     *W = linesearch->vec_sol_new;
14936a388c36SPeter Brune   }
14946a388c36SPeter Brune   if (G) {
14954f572ea9SToby Isaac     PetscAssertPointer(G, 6);
14966a388c36SPeter Brune     *G = linesearch->vec_func_new;
14976a388c36SPeter Brune   }
14983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14996a388c36SPeter Brune }
15006a388c36SPeter Brune 
1501f40b411bSPeter Brune /*@
1502f6dfbefdSBarry Smith   SNESLineSearchSetVecs - Sets the vectors on the `SNESLineSearch` context
1503f6dfbefdSBarry Smith 
1504c3339decSBarry Smith   Logically Collective
1505f40b411bSPeter Brune 
1506f40b411bSPeter Brune   Input Parameters:
1507420bcc1bSBarry Smith + linesearch - the line search context
15086232e825SPeter Brune . X          - Solution vector
15096232e825SPeter Brune . F          - Function vector
15106232e825SPeter Brune . Y          - Search direction vector
15116232e825SPeter Brune . W          - Solution work vector
15126232e825SPeter Brune - G          - Function work vector
1513f40b411bSPeter Brune 
1514420bcc1bSBarry Smith   Level: developer
1515f40b411bSPeter Brune 
1516420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetNorms()`, `SNESLineSearchGetVecs()`
1517f40b411bSPeter Brune @*/
1518d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch, Vec X, Vec F, Vec Y, Vec W, Vec G)
1519d71ae5a4SJacob Faibussowitsch {
15206a388c36SPeter Brune   PetscFunctionBegin;
1521f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
15226a388c36SPeter Brune   if (X) {
15236a388c36SPeter Brune     PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
15246a388c36SPeter Brune     linesearch->vec_sol = X;
15256a388c36SPeter Brune   }
15266a388c36SPeter Brune   if (F) {
15276a388c36SPeter Brune     PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
15286a388c36SPeter Brune     linesearch->vec_func = F;
15296a388c36SPeter Brune   }
15306a388c36SPeter Brune   if (Y) {
15316a388c36SPeter Brune     PetscValidHeaderSpecific(Y, VEC_CLASSID, 4);
15326a388c36SPeter Brune     linesearch->vec_update = Y;
15336a388c36SPeter Brune   }
15346a388c36SPeter Brune   if (W) {
15356a388c36SPeter Brune     PetscValidHeaderSpecific(W, VEC_CLASSID, 5);
15366a388c36SPeter Brune     linesearch->vec_sol_new = W;
15376a388c36SPeter Brune   }
15386a388c36SPeter Brune   if (G) {
15396a388c36SPeter Brune     PetscValidHeaderSpecific(G, VEC_CLASSID, 6);
15406a388c36SPeter Brune     linesearch->vec_func_new = G;
15416a388c36SPeter Brune   }
15423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15436a388c36SPeter Brune }
15446a388c36SPeter Brune 
1545cc4c1da9SBarry Smith /*@
1546f1c6b773SPeter Brune   SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1547f6dfbefdSBarry Smith   `SNESLineSearch` options in the database.
1548e7058c64SPeter Brune 
1549c3339decSBarry Smith   Logically Collective
1550e7058c64SPeter Brune 
1551e7058c64SPeter Brune   Input Parameters:
1552f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context
1553e7058c64SPeter Brune - prefix     - the prefix to prepend to all option names
1554e7058c64SPeter Brune 
155520f4b53cSBarry Smith   Level: advanced
155620f4b53cSBarry Smith 
1557f6dfbefdSBarry Smith   Note:
1558e7058c64SPeter Brune   A hyphen (-) must NOT be given at the beginning of the prefix name.
1559e7058c64SPeter Brune   The first character of all runtime options is AUTOMATICALLY the hyphen.
1560e7058c64SPeter Brune 
1561420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch()`, `SNESLineSearchSetFromOptions()`, `SNESGetOptionsPrefix()`
1562e7058c64SPeter Brune @*/
1563d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch, const char prefix[])
1564d71ae5a4SJacob Faibussowitsch {
1565e7058c64SPeter Brune   PetscFunctionBegin;
1566f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
15679566063dSJacob Faibussowitsch   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)linesearch, prefix));
15683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1569e7058c64SPeter Brune }
1570e7058c64SPeter Brune 
1571cc4c1da9SBarry Smith /*@
1572f6dfbefdSBarry Smith   SNESLineSearchGetOptionsPrefix - Gets the prefix used for searching for all
1573f1c6b773SPeter Brune   SNESLineSearch options in the database.
1574e7058c64SPeter Brune 
1575e7058c64SPeter Brune   Not Collective
1576e7058c64SPeter Brune 
1577e7058c64SPeter Brune   Input Parameter:
1578f6dfbefdSBarry Smith . linesearch - the `SNESLineSearch` context
1579e7058c64SPeter Brune 
1580e7058c64SPeter Brune   Output Parameter:
1581e7058c64SPeter Brune . prefix - pointer to the prefix string used
1582e7058c64SPeter Brune 
1583e7058c64SPeter Brune   Level: advanced
1584e7058c64SPeter Brune 
1585e4094ef1SJacob Faibussowitsch   Fortran Notes:
158620f4b53cSBarry Smith   The user should pass in a string 'prefix' of
158720f4b53cSBarry Smith   sufficient length to hold the prefix.
158820f4b53cSBarry Smith 
1589420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESAppendOptionsPrefix()`
1590e7058c64SPeter Brune @*/
1591d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch, const char *prefix[])
1592d71ae5a4SJacob Faibussowitsch {
1593e7058c64SPeter Brune   PetscFunctionBegin;
1594f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
15959566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)linesearch, prefix));
15963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1597e7058c64SPeter Brune }
1598bf7f4e0aSPeter Brune 
15998d359177SBarry Smith /*@C
1600f6dfbefdSBarry Smith   SNESLineSearchSetWorkVecs - Sets work vectors for the line search.
1601f40b411bSPeter Brune 
1602d8d19677SJose E. Roman   Input Parameters:
1603f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context
1604f40b411bSPeter Brune - nwork      - the number of work vectors
1605f40b411bSPeter Brune 
1606f40b411bSPeter Brune   Level: developer
1607f40b411bSPeter Brune 
1608420bcc1bSBarry Smith   Developer Note:
1609420bcc1bSBarry Smith   This is called from within the set up routines for each of the line search types `SNESLineSearchType`
1610420bcc1bSBarry Smith 
1611420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESSetWorkVecs()`
1612f40b411bSPeter Brune @*/
1613d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork)
1614d71ae5a4SJacob Faibussowitsch {
1615bf7f4e0aSPeter Brune   PetscFunctionBegin;
16160fdf79fbSJacob Faibussowitsch   PetscCheck(linesearch->vec_sol, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
16179566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work));
16183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1619bf7f4e0aSPeter Brune }
1620bf7f4e0aSPeter Brune 
1621f40b411bSPeter Brune /*@
1622422a814eSBarry Smith   SNESLineSearchGetReason - Gets the success/failure status of the last line search application
1623f40b411bSPeter Brune 
16242fe279fdSBarry Smith   Input Parameter:
1625420bcc1bSBarry Smith . linesearch - the line search context
1626f40b411bSPeter Brune 
16272fe279fdSBarry Smith   Output Parameter:
1628422a814eSBarry Smith . result - The success or failure status
1629f40b411bSPeter Brune 
163020f4b53cSBarry Smith   Level: developer
163120f4b53cSBarry Smith 
1632f6dfbefdSBarry Smith   Note:
1633420bcc1bSBarry Smith   This is typically called after `SNESLineSearchApply()` in order to determine if the line search failed
1634420bcc1bSBarry Smith   (and set into the `SNES` convergence accordingly).
1635cd7522eaSPeter Brune 
1636420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetReason()`, `SNESLineSearchReason`
1637f40b411bSPeter Brune @*/
1638d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result)
1639d71ae5a4SJacob Faibussowitsch {
1640bf7f4e0aSPeter Brune   PetscFunctionBegin;
1641f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
16424f572ea9SToby Isaac   PetscAssertPointer(result, 2);
1643422a814eSBarry Smith   *result = linesearch->result;
16443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1645bf7f4e0aSPeter Brune }
1646bf7f4e0aSPeter Brune 
1647f40b411bSPeter Brune /*@
1648420bcc1bSBarry Smith   SNESLineSearchSetReason - Sets the success/failure status of the line search application
1649f40b411bSPeter Brune 
16505d83a8b1SBarry Smith   Logically Collective; No Fortran Support
16515d83a8b1SBarry Smith 
1652f40b411bSPeter Brune   Input Parameters:
1653420bcc1bSBarry Smith + linesearch - the line search context
1654422a814eSBarry Smith - result     - The success or failure status
1655f40b411bSPeter Brune 
165620f4b53cSBarry Smith   Level: developer
165720f4b53cSBarry Smith 
1658f6dfbefdSBarry Smith   Note:
1659420bcc1bSBarry Smith   This is typically called in a `SNESLineSearchType` implementation of `SNESLineSearchApply()` or a `SNESLINESEARCHSHELL` implementation to set
1660cd7522eaSPeter Brune   the success or failure of the line search method.
1661cd7522eaSPeter Brune 
1662420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchReason`, `SNESLineSearchGetSResult()`
1663f40b411bSPeter Brune @*/
1664d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result)
1665d71ae5a4SJacob Faibussowitsch {
16666a388c36SPeter Brune   PetscFunctionBegin;
16675fd66863SKarl Rupp   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1668422a814eSBarry Smith   linesearch->result = result;
16693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16706a388c36SPeter Brune }
16716a388c36SPeter Brune 
1672ceaaa498SBarry Smith // PetscClangLinter pragma disable: -fdoc-param-list-func-parameter-documentation
16739bd66eb0SPeter Brune /*@C
1674f1c6b773SPeter Brune   SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.
16759bd66eb0SPeter Brune 
1676c3339decSBarry Smith   Logically Collective
1677f6dfbefdSBarry Smith 
16789bd66eb0SPeter Brune   Input Parameters:
1679ceaaa498SBarry Smith + linesearch  - the linesearch object
16808434afd1SBarry Smith . projectfunc - function for projecting the function to the bounds, see `SNESLineSearchVIProjectFn` for calling sequence
16818434afd1SBarry Smith - normfunc    - function for computing the norm of an active set, see `SNESLineSearchVINormFn` for calling sequence
16829bd66eb0SPeter Brune 
1683f6dfbefdSBarry Smith   Level: advanced
16849bd66eb0SPeter Brune 
1685ceaaa498SBarry Smith   Notes:
168620f4b53cSBarry Smith   The VI solvers require projection of the solution to the feasible set.  `projectfunc` should implement this.
168720f4b53cSBarry Smith 
168820f4b53cSBarry Smith   The VI solvers require special evaluation of the function norm such that the norm is only calculated
168920f4b53cSBarry Smith   on the inactive set.  This should be implemented by `normfunc`.
169020f4b53cSBarry Smith 
16919bcc50f1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetVIFunctions()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()`,
16928434afd1SBarry Smith           `SNESLineSearchVIProjectFn`, `SNESLineSearchVINormFn`
16939bd66eb0SPeter Brune @*/
16948434afd1SBarry Smith PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFn *projectfunc, SNESLineSearchVINormFn *normfunc)
1695d71ae5a4SJacob Faibussowitsch {
16969bd66eb0SPeter Brune   PetscFunctionBegin;
1697f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
16989bd66eb0SPeter Brune   if (projectfunc) linesearch->ops->viproject = projectfunc;
16999bd66eb0SPeter Brune   if (normfunc) linesearch->ops->vinorm = normfunc;
17003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17019bd66eb0SPeter Brune }
17029bd66eb0SPeter Brune 
17039bd66eb0SPeter Brune /*@C
1704f1c6b773SPeter Brune   SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.
17059bd66eb0SPeter Brune 
1706f6dfbefdSBarry Smith   Not Collective
1707f6dfbefdSBarry Smith 
1708f899ff85SJose E. Roman   Input Parameter:
1709f6dfbefdSBarry Smith . linesearch - the line search context, obtain with `SNESGetLineSearch()`
17109bd66eb0SPeter Brune 
17119bd66eb0SPeter Brune   Output Parameters:
17128434afd1SBarry Smith + projectfunc - function for projecting the function to the bounds, see `SNESLineSearchVIProjectFn` for calling sequence
17138434afd1SBarry Smith - normfunc    - function for computing the norm of an active set, see `SNESLineSearchVINormFn ` for calling sequence
17149bd66eb0SPeter Brune 
1715f6dfbefdSBarry Smith   Level: advanced
17169bd66eb0SPeter Brune 
17179bcc50f1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetVIFunctions()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()`,
17188434afd1SBarry Smith           `SNESLineSearchVIProjectFn`, `SNESLineSearchVINormFn`
17199bd66eb0SPeter Brune @*/
17208434afd1SBarry Smith PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFn **projectfunc, SNESLineSearchVINormFn **normfunc)
1721d71ae5a4SJacob Faibussowitsch {
17229bd66eb0SPeter Brune   PetscFunctionBegin;
17239bd66eb0SPeter Brune   if (projectfunc) *projectfunc = linesearch->ops->viproject;
17249bd66eb0SPeter Brune   if (normfunc) *normfunc = linesearch->ops->vinorm;
17253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17269bd66eb0SPeter Brune }
17279bd66eb0SPeter Brune 
1728bf7f4e0aSPeter Brune /*@C
1729420bcc1bSBarry Smith   SNESLineSearchRegister - register a line search type `SNESLineSearchType`
1730ceaaa498SBarry Smith 
1731cc4c1da9SBarry Smith   Logically Collective, No Fortran Support
1732cc4c1da9SBarry Smith 
1733ceaaa498SBarry Smith   Input Parameters:
1734420bcc1bSBarry Smith + sname    - name of the `SNESLineSearchType()`
1735ceaaa498SBarry Smith - function - the creation function for that type
1736ceaaa498SBarry Smith 
1737ceaaa498SBarry Smith   Calling sequence of `function`:
1738ceaaa498SBarry Smith . ls - the line search context
1739bf7f4e0aSPeter Brune 
1740bf7f4e0aSPeter Brune   Level: advanced
1741f6dfbefdSBarry Smith 
1742420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchSetType()`
1743bf7f4e0aSPeter Brune @*/
1744ceaaa498SBarry Smith PetscErrorCode SNESLineSearchRegister(const char sname[], PetscErrorCode (*function)(SNESLineSearch ls))
1745d71ae5a4SJacob Faibussowitsch {
1746bf7f4e0aSPeter Brune   PetscFunctionBegin;
17479566063dSJacob Faibussowitsch   PetscCall(SNESInitializePackage());
17489566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&SNESLineSearchList, sname, function));
17493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1750bf7f4e0aSPeter Brune }
1751