xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision f6dfbefd03961ab3be6f06be75c96cbf27a49667)
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 /*@
10*f6dfbefdSBarry Smith    SNESLineSearchMonitorCancel - Clears all the monitor functions for a `SNESLineSearch` object.
11dcf2fd19SBarry Smith 
12*f6dfbefdSBarry Smith    Logically Collective on ls
13dcf2fd19SBarry Smith 
14dcf2fd19SBarry Smith    Input Parameters:
15*f6dfbefdSBarry 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
19*f6dfbefdSBarry Smith     into a code by calls to `SNESLineSearchMonitorSet()`, but does not cancel those
20dcf2fd19SBarry Smith     set via the options database
21dcf2fd19SBarry Smith 
22dcf2fd19SBarry Smith    Notes:
23*f6dfbefdSBarry Smith    There is no way to clear one specific monitor from a `SNESLineSearch` object.
24dcf2fd19SBarry Smith 
25*f6dfbefdSBarry Smith    This does not clear the monitor set with `SNESLineSearchSetDefaultMonitor()` use `SNESLineSearchSetDefaultMonitor`(ls,NULL) to cancel
26dcf2fd19SBarry Smith    that one.
27dcf2fd19SBarry Smith 
28*f6dfbefdSBarry Smith    Level: advanced
29dcf2fd19SBarry Smith 
30*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorDefault()`, `SNESLineSearchMonitorSet()`
31dcf2fd19SBarry Smith @*/
329371c9d4SSatish Balay PetscErrorCode SNESLineSearchMonitorCancel(SNESLineSearch ls) {
33dcf2fd19SBarry Smith   PetscInt i;
34dcf2fd19SBarry Smith 
35dcf2fd19SBarry Smith   PetscFunctionBegin;
36dcf2fd19SBarry Smith   PetscValidHeaderSpecific(ls, SNESLINESEARCH_CLASSID, 1);
37dcf2fd19SBarry Smith   for (i = 0; i < ls->numbermonitors; i++) {
3848a46eb9SPierre Jolivet     if (ls->monitordestroy[i]) PetscCall((*ls->monitordestroy[i])(&ls->monitorcontext[i]));
39dcf2fd19SBarry Smith   }
40dcf2fd19SBarry Smith   ls->numbermonitors = 0;
41dcf2fd19SBarry Smith   PetscFunctionReturn(0);
42dcf2fd19SBarry Smith }
43dcf2fd19SBarry Smith 
44dcf2fd19SBarry Smith /*@
45dcf2fd19SBarry Smith    SNESLineSearchMonitor - runs the user provided monitor routines, if they exist
46dcf2fd19SBarry Smith 
47*f6dfbefdSBarry Smith    Collective on ls
48dcf2fd19SBarry Smith 
49dcf2fd19SBarry Smith    Input Parameters:
50dcf2fd19SBarry Smith .  ls - the linesearch object
51dcf2fd19SBarry Smith 
52*f6dfbefdSBarry Smith    Note:
53dcf2fd19SBarry Smith    This routine is called by the SNES implementations.
54dcf2fd19SBarry Smith    It does not typically need to be called by the user.
55dcf2fd19SBarry Smith 
56dcf2fd19SBarry Smith    Level: developer
57dcf2fd19SBarry Smith 
58*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorSet()`
59dcf2fd19SBarry Smith @*/
609371c9d4SSatish Balay PetscErrorCode SNESLineSearchMonitor(SNESLineSearch ls) {
61dcf2fd19SBarry Smith   PetscInt i, n = ls->numbermonitors;
62dcf2fd19SBarry Smith 
63dcf2fd19SBarry Smith   PetscFunctionBegin;
6448a46eb9SPierre Jolivet   for (i = 0; i < n; i++) PetscCall((*ls->monitorftns[i])(ls, ls->monitorcontext[i]));
65dcf2fd19SBarry Smith   PetscFunctionReturn(0);
66dcf2fd19SBarry Smith }
67dcf2fd19SBarry Smith 
68dcf2fd19SBarry Smith /*@C
69dcf2fd19SBarry Smith    SNESLineSearchMonitorSet - Sets an ADDITIONAL function that is to be used at every
70dcf2fd19SBarry Smith    iteration of the nonlinear solver to display the iteration's
71dcf2fd19SBarry Smith    progress.
72dcf2fd19SBarry Smith 
73*f6dfbefdSBarry Smith    Logically Collective on ls
74dcf2fd19SBarry Smith 
75dcf2fd19SBarry Smith    Input Parameters:
76*f6dfbefdSBarry Smith +  ls - the `SNESLineSearch` context
77dcf2fd19SBarry Smith .  f - the monitor function
78dcf2fd19SBarry Smith .  mctx - [optional] user-defined context for private data for the
79dcf2fd19SBarry Smith           monitor routine (use NULL if no context is desired)
80dcf2fd19SBarry Smith -  monitordestroy - [optional] routine that frees monitor context
81dcf2fd19SBarry Smith           (may be NULL)
82dcf2fd19SBarry Smith 
83*f6dfbefdSBarry Smith    Note:
84dcf2fd19SBarry Smith    Several different monitoring routines may be set by calling
85*f6dfbefdSBarry Smith    `SNESLineSearchMonitorSet()` multiple times; all will be called in the
86dcf2fd19SBarry Smith    order in which they were set.
87dcf2fd19SBarry Smith 
88*f6dfbefdSBarry Smith    Fortran Note:
89*f6dfbefdSBarry Smith    Only a single monitor function can be set for each `SNESLineSearch` object
90dcf2fd19SBarry Smith 
91dcf2fd19SBarry Smith    Level: intermediate
92dcf2fd19SBarry Smith 
93*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorDefault()`, `SNESLineSearchMonitorCancel()`
94dcf2fd19SBarry Smith @*/
959371c9d4SSatish Balay PetscErrorCode SNESLineSearchMonitorSet(SNESLineSearch ls, PetscErrorCode (*f)(SNESLineSearch, void *), void *mctx, PetscErrorCode (*monitordestroy)(void **)) {
9678064530SBarry Smith   PetscInt  i;
9778064530SBarry Smith   PetscBool identical;
9878064530SBarry Smith 
99dcf2fd19SBarry Smith   PetscFunctionBegin;
100dcf2fd19SBarry Smith   PetscValidHeaderSpecific(ls, SNESLINESEARCH_CLASSID, 1);
10178064530SBarry Smith   for (i = 0; i < ls->numbermonitors; i++) {
1029566063dSJacob Faibussowitsch     PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))f, mctx, monitordestroy, (PetscErrorCode(*)(void))ls->monitorftns[i], ls->monitorcontext[i], ls->monitordestroy[i], &identical));
10378064530SBarry Smith     if (identical) PetscFunctionReturn(0);
10478064530SBarry Smith   }
1055f80ce2aSJacob Faibussowitsch   PetscCheck(ls->numbermonitors < MAXSNESLSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many monitors set");
106dcf2fd19SBarry Smith   ls->monitorftns[ls->numbermonitors]      = f;
107dcf2fd19SBarry Smith   ls->monitordestroy[ls->numbermonitors]   = monitordestroy;
108dcf2fd19SBarry Smith   ls->monitorcontext[ls->numbermonitors++] = (void *)mctx;
109dcf2fd19SBarry Smith   PetscFunctionReturn(0);
110dcf2fd19SBarry Smith }
111dcf2fd19SBarry Smith 
112dcf2fd19SBarry Smith /*@C
113*f6dfbefdSBarry Smith    SNESLineSearchMonitorSolutionUpdate - Monitors each update of the function value the linesearch tries
114dcf2fd19SBarry Smith 
115*f6dfbefdSBarry Smith    Collective on ls
116dcf2fd19SBarry Smith 
117dcf2fd19SBarry Smith    Input Parameters:
118*f6dfbefdSBarry Smith +  ls - the `SNES` linesearch object
119*f6dfbefdSBarry Smith -  vf - the context for the monitor, in this case it is an `PetscViewerAndFormat`
120dcf2fd19SBarry Smith 
121dcf2fd19SBarry Smith    Level: intermediate
122dcf2fd19SBarry Smith 
123*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESMonitorSet()`, `SNESMonitorSolution()`
124dcf2fd19SBarry Smith @*/
1259371c9d4SSatish Balay PetscErrorCode SNESLineSearchMonitorSolutionUpdate(SNESLineSearch ls, PetscViewerAndFormat *vf) {
126d12e167eSBarry Smith   PetscViewer viewer = vf->viewer;
127dcf2fd19SBarry Smith   Vec         Y, W, G;
128dcf2fd19SBarry Smith 
129dcf2fd19SBarry Smith   PetscFunctionBegin;
1309566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetVecs(ls, NULL, NULL, &Y, &W, &G));
1319566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(viewer, vf->format));
1329566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted update to solution \n"));
1339566063dSJacob Faibussowitsch   PetscCall(VecView(Y, viewer));
1349566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted new solution \n"));
1359566063dSJacob Faibussowitsch   PetscCall(VecView(W, viewer));
1369566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted updated function value\n"));
1379566063dSJacob Faibussowitsch   PetscCall(VecView(G, viewer));
1389566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(viewer));
139dcf2fd19SBarry Smith   PetscFunctionReturn(0);
140dcf2fd19SBarry Smith }
141dcf2fd19SBarry Smith 
142f40b411bSPeter Brune /*@
143cd7522eaSPeter Brune    SNESLineSearchCreate - Creates the line search context.
144f40b411bSPeter Brune 
145*f6dfbefdSBarry Smith    Logically Collective
146f40b411bSPeter Brune 
147f40b411bSPeter Brune    Input Parameters:
148*f6dfbefdSBarry Smith .  comm - MPI communicator for the line search (typically from the associated `SNES` context).
149f40b411bSPeter Brune 
150f40b411bSPeter Brune    Output Parameters:
1518e557f58SPeter Brune .  outlinesearch - the new linesearch context
152f40b411bSPeter Brune 
153162e0bf5SPeter Brune    Level: developer
154162e0bf5SPeter Brune 
155*f6dfbefdSBarry Smith    Note:
156*f6dfbefdSBarry Smith    The preferred calling sequence for users is to use `SNESGetLineSearch()` to acquire the `SNESLineSearch` instance
157*f6dfbefdSBarry Smith    already associated with the `SNES`.
158f40b411bSPeter Brune 
159*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `LineSearchDestroy()`, `SNESGetLineSearch()`
160f40b411bSPeter Brune @*/
1619371c9d4SSatish Balay PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch) {
162f1c6b773SPeter Brune   SNESLineSearch linesearch;
163bf388a1fSBarry Smith 
164bf7f4e0aSPeter Brune   PetscFunctionBegin;
165ea5d4fccSPeter Brune   PetscValidPointer(outlinesearch, 2);
1669566063dSJacob Faibussowitsch   PetscCall(SNESInitializePackage());
1670298fd71SBarry Smith   *outlinesearch = NULL;
168f5af7f23SKarl Rupp 
1699566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(linesearch, SNESLINESEARCH_CLASSID, "SNESLineSearch", "Linesearch", "SNESLineSearch", comm, SNESLineSearchDestroy, SNESLineSearchView));
170bf7f4e0aSPeter Brune 
1710298fd71SBarry Smith   linesearch->vec_sol_new  = NULL;
1720298fd71SBarry Smith   linesearch->vec_func_new = NULL;
1730298fd71SBarry Smith   linesearch->vec_sol      = NULL;
1740298fd71SBarry Smith   linesearch->vec_func     = NULL;
1750298fd71SBarry Smith   linesearch->vec_update   = NULL;
1769bd66eb0SPeter Brune 
177bf7f4e0aSPeter Brune   linesearch->lambda       = 1.0;
178bf7f4e0aSPeter Brune   linesearch->fnorm        = 1.0;
179bf7f4e0aSPeter Brune   linesearch->ynorm        = 1.0;
180bf7f4e0aSPeter Brune   linesearch->xnorm        = 1.0;
181422a814eSBarry Smith   linesearch->result       = SNES_LINESEARCH_SUCCEEDED;
182bf7f4e0aSPeter Brune   linesearch->norms        = PETSC_TRUE;
183bf7f4e0aSPeter Brune   linesearch->keeplambda   = PETSC_FALSE;
184bf7f4e0aSPeter Brune   linesearch->damping      = 1.0;
185bf7f4e0aSPeter Brune   linesearch->maxstep      = 1e8;
186bf7f4e0aSPeter Brune   linesearch->steptol      = 1e-12;
187516fe3c3SPeter Brune   linesearch->rtol         = 1e-8;
188516fe3c3SPeter Brune   linesearch->atol         = 1e-15;
189516fe3c3SPeter Brune   linesearch->ltol         = 1e-8;
1900298fd71SBarry Smith   linesearch->precheckctx  = NULL;
1910298fd71SBarry Smith   linesearch->postcheckctx = NULL;
19259405d5eSPeter Brune   linesearch->max_its      = 1;
193bf7f4e0aSPeter Brune   linesearch->setupcalled  = PETSC_FALSE;
1943add74b1SBarry Smith   linesearch->monitor      = NULL;
195bf7f4e0aSPeter Brune   *outlinesearch           = linesearch;
196bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
197bf7f4e0aSPeter Brune }
198bf7f4e0aSPeter Brune 
199f40b411bSPeter Brune /*@
20078bcb3b5SPeter Brune    SNESLineSearchSetUp - Prepares the line search for being applied by allocating
20178bcb3b5SPeter Brune    any required vectors.
202f40b411bSPeter Brune 
203*f6dfbefdSBarry Smith    Collective on linesearch
204f40b411bSPeter Brune 
205f40b411bSPeter Brune    Input Parameters:
206*f6dfbefdSBarry Smith .  linesearch - The `SNESLineSearch` instance.
207f40b411bSPeter Brune 
208*f6dfbefdSBarry Smith    Note:
209*f6dfbefdSBarry Smith    For most cases, this needn't be called by users or outside of `SNESLineSearchApply()`.
210cd7522eaSPeter Brune    The only current case where this is called outside of this is for the VI
21178bcb3b5SPeter Brune    solvers, which modify the solution and work vectors before the first call
212*f6dfbefdSBarry Smith    of `SNESLineSearchApply()`, requiring the `SNESLineSearch` work vectors to be
213cd7522eaSPeter Brune    allocated upfront.
214cd7522eaSPeter Brune 
21578bcb3b5SPeter Brune    Level: advanced
216f40b411bSPeter Brune 
217*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchReset()`
218f40b411bSPeter Brune @*/
219f40b411bSPeter Brune 
2209371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch) {
221bf7f4e0aSPeter Brune   PetscFunctionBegin;
22248a46eb9SPierre Jolivet   if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC));
223bf7f4e0aSPeter Brune   if (!linesearch->setupcalled) {
22448a46eb9SPierre Jolivet     if (!linesearch->vec_sol_new) PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new));
22548a46eb9SPierre Jolivet     if (!linesearch->vec_func_new) PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new));
226dbbe0bcdSBarry Smith     if (linesearch->ops->setup) PetscUseTypeMethod(linesearch, setup);
2279566063dSJacob Faibussowitsch     if (!linesearch->ops->snesfunc) PetscCall(SNESLineSearchSetFunction(linesearch, SNESComputeFunction));
228bf7f4e0aSPeter Brune     linesearch->lambda      = linesearch->damping;
229bf7f4e0aSPeter Brune     linesearch->setupcalled = PETSC_TRUE;
230bf7f4e0aSPeter Brune   }
231bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
232bf7f4e0aSPeter Brune }
233bf7f4e0aSPeter Brune 
234f40b411bSPeter Brune /*@
235*f6dfbefdSBarry Smith    SNESLineSearchReset - Undoes the `SNESLineSearchSetUp()` and deletes any `Vec`s or `Mat`s allocated by the line search.
236f40b411bSPeter Brune 
237*f6dfbefdSBarry Smith    Collective on linesearch
238f40b411bSPeter Brune 
239f40b411bSPeter Brune    Input Parameters:
240*f6dfbefdSBarry Smith .  linesearch - The `SNESLineSearch` instance.
241f40b411bSPeter Brune 
242*f6dfbefdSBarry Smith    Note:
243*f6dfbefdSBarry Smith     Usually only called by `SNESReset()`
244f190f2fcSBarry Smith 
245f190f2fcSBarry Smith    Level: developer
246f40b411bSPeter Brune 
247*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchSetUp()`
248f40b411bSPeter Brune @*/
249f40b411bSPeter Brune 
2509371c9d4SSatish Balay PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch) {
251bf7f4e0aSPeter Brune   PetscFunctionBegin;
252dbbe0bcdSBarry Smith   if (linesearch->ops->reset) PetscUseTypeMethod(linesearch, reset);
253f5af7f23SKarl Rupp 
2549566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&linesearch->vec_sol_new));
2559566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&linesearch->vec_func_new));
256bf7f4e0aSPeter Brune 
2579566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(linesearch->nwork, &linesearch->work));
258f5af7f23SKarl Rupp 
259bf7f4e0aSPeter Brune   linesearch->nwork       = 0;
260bf7f4e0aSPeter Brune   linesearch->setupcalled = PETSC_FALSE;
261bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
262bf7f4e0aSPeter Brune }
263bf7f4e0aSPeter Brune 
264ed07d7d7SPeter Brune /*@C
265*f6dfbefdSBarry Smith    SNESLineSearchSetFunction - Sets the function evaluation used by the `SNES` line search
266*f6dfbefdSBarry Smith `
267ed07d7d7SPeter Brune    Input Parameters:
268*f6dfbefdSBarry Smith .  linesearch - the `SNESLineSearch` context
269*f6dfbefdSBarry Smith +  func       - function evaluation routine, this is usually the function provided with `SNESSetFunction()`
270ed07d7d7SPeter Brune 
271ed07d7d7SPeter Brune    Level: developer
272ed07d7d7SPeter Brune 
273*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESSetFunction()`
274ed07d7d7SPeter Brune @*/
2759371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES, Vec, Vec)) {
276ed07d7d7SPeter Brune   PetscFunctionBegin;
277ed07d7d7SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
278ed07d7d7SPeter Brune   linesearch->ops->snesfunc = func;
279ed07d7d7SPeter Brune   PetscFunctionReturn(0);
280ed07d7d7SPeter Brune }
281ed07d7d7SPeter Brune 
28286d74e61SPeter Brune /*@C
283f190f2fcSBarry Smith    SNESLineSearchSetPreCheck - Sets a user function that is called after the initial search direction has been computed but
284df3898eeSBarry Smith          before the line search routine has been applied. Allows the user to adjust the result of (usually a linear solve) that
285f190f2fcSBarry Smith          determined the search direction.
28686d74e61SPeter Brune 
287*f6dfbefdSBarry Smith    Logically Collective on linesearch
28886d74e61SPeter Brune 
28986d74e61SPeter Brune    Input Parameters:
290*f6dfbefdSBarry Smith +  linesearch - the `SNESLineSearch` context
291*f6dfbefdSBarry Smith .  func - [optional] function evaluation routine, see `SNESLineSearchPreCheck()` for the calling sequence
292f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
29386d74e61SPeter Brune 
29486d74e61SPeter Brune    Level: intermediate
29586d74e61SPeter Brune 
296*f6dfbefdSBarry Smith    Note:
297f0b84518SBarry Smith    Use `SNESLineSearchSetPostCheck()` to change the step after the line search.
298f0b84518SBarry Smith    search is complete.
299f0b84518SBarry Smith 
300f0b84518SBarry Smith    Use `SNESVISetVariableBounds()` and `SNESVISetComputeVariableBounds()` to cause `SNES` to automatically control the ranges of variables allowed.
301f0b84518SBarry Smith 
302f0b84518SBarry Smith .seealso: `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()`,
303f0b84518SBarry Smith           `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()
304f0b84518SBarry Smith 
30586d74e61SPeter Brune @*/
3069371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch, Vec, Vec, PetscBool *, void *), void *ctx) {
3079bd66eb0SPeter Brune   PetscFunctionBegin;
308f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
309f190f2fcSBarry Smith   if (func) linesearch->ops->precheck = func;
31086d74e61SPeter Brune   if (ctx) linesearch->precheckctx = ctx;
31186d74e61SPeter Brune   PetscFunctionReturn(0);
31286d74e61SPeter Brune }
31386d74e61SPeter Brune 
31486d74e61SPeter Brune /*@C
31553deb899SBarry Smith    SNESLineSearchGetPreCheck - Gets the pre-check function for the line search routine.
31686d74e61SPeter Brune 
317f899ff85SJose E. Roman    Input Parameter:
318*f6dfbefdSBarry Smith .  linesearch - the `SNESLineSearch` context
31986d74e61SPeter Brune 
32086d74e61SPeter Brune    Output Parameters:
321*f6dfbefdSBarry Smith +  func       - [optional] function evaluation routine, see `SNESLineSearchPreCheck` for calling sequence
322f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
32386d74e61SPeter Brune 
32486d74e61SPeter Brune    Level: intermediate
32586d74e61SPeter Brune 
326*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchSetPostCheck()`
32786d74e61SPeter Brune @*/
3289371c9d4SSatish Balay PetscErrorCode SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch, Vec, Vec, PetscBool *, void *), void **ctx) {
3299bd66eb0SPeter Brune   PetscFunctionBegin;
330f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
331f190f2fcSBarry Smith   if (func) *func = linesearch->ops->precheck;
33286d74e61SPeter Brune   if (ctx) *ctx = linesearch->precheckctx;
33386d74e61SPeter Brune   PetscFunctionReturn(0);
33486d74e61SPeter Brune }
33586d74e61SPeter Brune 
33686d74e61SPeter Brune /*@C
337f190f2fcSBarry Smith    SNESLineSearchSetPostCheck - Sets a user function that is called after the line search has been applied to determine the step
338f190f2fcSBarry Smith        direction and length. Allows the user a chance to change or override the decision of the line search routine
33986d74e61SPeter Brune 
340*f6dfbefdSBarry Smith    Logically Collective linesearch
34186d74e61SPeter Brune 
34286d74e61SPeter Brune    Input Parameters:
343*f6dfbefdSBarry Smith +  linesearch - the `SNESLineSearch` context
344*f6dfbefdSBarry Smith .  func - [optional] function evaluation routine, see `SNESLineSearchPostCheck`  for the calling sequence
345f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
34686d74e61SPeter Brune 
34786d74e61SPeter Brune    Level: intermediate
34886d74e61SPeter Brune 
349f0b84518SBarry Smith    Notes:
350f0b84518SBarry Smith    Use `SNESLineSearchSetPreCheck()` to change the step before the line search.
351f0b84518SBarry Smith    search is complete.
352f0b84518SBarry Smith 
353f0b84518SBarry Smith    Use `SNESVISetVariableBounds()` and `SNESVISetComputeVariableBounds()` to cause `SNES` to automatically control the ranges of variables allowed.
354f0b84518SBarry Smith 
355*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchGetPostCheck()`,
356f0b84518SBarry Smith           `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()
35786d74e61SPeter Brune @*/
3589371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx) {
35986d74e61SPeter Brune   PetscFunctionBegin;
360f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
361f190f2fcSBarry Smith   if (func) linesearch->ops->postcheck = func;
36286d74e61SPeter Brune   if (ctx) linesearch->postcheckctx = ctx;
36386d74e61SPeter Brune   PetscFunctionReturn(0);
36486d74e61SPeter Brune }
36586d74e61SPeter Brune 
36686d74e61SPeter Brune /*@C
367f1c6b773SPeter Brune    SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine.
36886d74e61SPeter Brune 
369f899ff85SJose E. Roman    Input Parameter:
370*f6dfbefdSBarry Smith .  linesearch - the `SNESLineSearch` context
37186d74e61SPeter Brune 
37286d74e61SPeter Brune    Output Parameters:
373*f6dfbefdSBarry Smith +  func - [optional] function evaluation routine, see for the calling sequence `SNESLineSearchPostCheck`
374f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
37586d74e61SPeter Brune 
37686d74e61SPeter Brune    Level: intermediate
37786d74e61SPeter Brune 
378*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`
37986d74e61SPeter Brune @*/
3809371c9d4SSatish Balay PetscErrorCode SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx) {
3819bd66eb0SPeter Brune   PetscFunctionBegin;
382f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
383f190f2fcSBarry Smith   if (func) *func = linesearch->ops->postcheck;
38486d74e61SPeter Brune   if (ctx) *ctx = linesearch->postcheckctx;
38586d74e61SPeter Brune   PetscFunctionReturn(0);
38686d74e61SPeter Brune }
38786d74e61SPeter Brune 
388f40b411bSPeter Brune /*@
389f1c6b773SPeter Brune    SNESLineSearchPreCheck - Prepares the line search for being applied.
390f40b411bSPeter Brune 
391*f6dfbefdSBarry Smith    Logically Collective on linesearch
392f40b411bSPeter Brune 
393f40b411bSPeter Brune    Input Parameters:
3947b1df9c1SPeter Brune +  linesearch - The linesearch instance.
3957b1df9c1SPeter Brune .  X - The current solution
3967b1df9c1SPeter Brune -  Y - The step direction
397f40b411bSPeter Brune 
398f40b411bSPeter Brune    Output Parameters:
3998e557f58SPeter Brune .  changed - Indicator that the precheck routine has changed anything
400f40b411bSPeter Brune 
401f0b84518SBarry Smith    Level: advanced
402f40b411bSPeter Brune 
403*f6dfbefdSBarry Smith    Note:
404*f6dfbefdSBarry Smith    This calls any function provided with `SNESLineSearchSetPreCheck()`
405*f6dfbefdSBarry Smith 
406*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`,  `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`,
407f0b84518SBarry Smith           `SNESLineSearchGetPostCheck()``
408f40b411bSPeter Brune @*/
4099371c9d4SSatish Balay PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed) {
410bf7f4e0aSPeter Brune   PetscFunctionBegin;
411bf7f4e0aSPeter Brune   *changed = PETSC_FALSE;
4126b2b7091SBarry Smith   if (linesearch->ops->precheck) {
413dbbe0bcdSBarry Smith     PetscUseTypeMethod(linesearch, precheck, X, Y, changed, linesearch->precheckctx);
41438bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch, *changed, 4);
415bf7f4e0aSPeter Brune   }
416bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
417bf7f4e0aSPeter Brune }
418bf7f4e0aSPeter Brune 
419f40b411bSPeter Brune /*@
420ef46b1a6SStefano Zampini    SNESLineSearchPostCheck - Hook to modify step direction or updated solution after a successful linesearch
421f40b411bSPeter Brune 
422*f6dfbefdSBarry Smith    Logically Collective on linesearch
423f40b411bSPeter Brune 
424f40b411bSPeter Brune    Input Parameters:
4257b1df9c1SPeter Brune +  linesearch - The linesearch context
4267b1df9c1SPeter Brune .  X - The last solution
4277b1df9c1SPeter Brune .  Y - The step direction
4287b1df9c1SPeter Brune -  W - The updated solution, W = X + lambda*Y for some lambda
429f40b411bSPeter Brune 
430f40b411bSPeter Brune    Output Parameters:
43178bcb3b5SPeter Brune +  changed_Y - Indicator if the direction Y has been changed.
43278bcb3b5SPeter Brune -  changed_W - Indicator if the new candidate solution W has been changed.
433f40b411bSPeter Brune 
434*f6dfbefdSBarry Smith    Note:
435*f6dfbefdSBarry Smith    This calls any function provided with `SNESLineSearchSetPreCheck()`
436*f6dfbefdSBarry Smith 
437f190f2fcSBarry Smith    Level: developer
438f40b411bSPeter Brune 
439db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPrecheck()`, `SNESLineSearchGetPrecheck()`
440f40b411bSPeter Brune @*/
4419371c9d4SSatish Balay PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W) {
442bf7f4e0aSPeter Brune   PetscFunctionBegin;
443bf7f4e0aSPeter Brune   *changed_Y = PETSC_FALSE;
444bf7f4e0aSPeter Brune   *changed_W = PETSC_FALSE;
4456b2b7091SBarry Smith   if (linesearch->ops->postcheck) {
446dbbe0bcdSBarry Smith     PetscUseTypeMethod(linesearch, postcheck, X, Y, W, changed_Y, changed_W, linesearch->postcheckctx);
44738bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch, *changed_Y, 5);
44838bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch, *changed_W, 6);
44986d74e61SPeter Brune   }
45086d74e61SPeter Brune   PetscFunctionReturn(0);
45186d74e61SPeter Brune }
45286d74e61SPeter Brune 
45386d74e61SPeter Brune /*@C
45486d74e61SPeter Brune    SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration
45586d74e61SPeter Brune 
456*f6dfbefdSBarry Smith    Logically Collective on linesearch
45786d74e61SPeter Brune 
4584165533cSJose E. Roman    Input Parameters:
45986d74e61SPeter Brune +  linesearch - linesearch context
46086d74e61SPeter Brune .  X - base state for this step
461907376e6SBarry Smith -  ctx - context for this function
46286d74e61SPeter Brune 
46397bb3fdcSJose E. Roman    Input/Output Parameter:
46497bb3fdcSJose E. Roman .  Y - correction, possibly modified
46597bb3fdcSJose E. Roman 
46697bb3fdcSJose E. Roman    Output Parameter:
46797bb3fdcSJose E. Roman .  changed - flag indicating that Y was modified
46886d74e61SPeter Brune 
46986d74e61SPeter Brune    Options Database Key:
470cd7522eaSPeter Brune +  -snes_linesearch_precheck_picard - activate this routine
471cd7522eaSPeter Brune -  -snes_linesearch_precheck_picard_angle - angle
47286d74e61SPeter Brune 
47386d74e61SPeter Brune    Level: advanced
47486d74e61SPeter Brune 
47586d74e61SPeter Brune    Notes:
476*f6dfbefdSBarry Smith    This function should be passed to `SNESLineSearchSetPreCheck()`
47786d74e61SPeter Brune 
47886d74e61SPeter Brune    The justification for this method involves the linear convergence of a Picard iteration
47986d74e61SPeter Brune    so the Picard linearization should be provided in place of the "Jacobian". This correction
48086d74e61SPeter Brune    is generally not useful when using a Newton linearization.
48186d74e61SPeter Brune 
48286d74e61SPeter Brune    Reference:
483*f6dfbefdSBarry Smith  . - * - Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology.
48486d74e61SPeter Brune 
485*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESSetPicard()`, `SNESGetLineSearch()`, `SNESLineSearchSetPreCheck()`
48686d74e61SPeter Brune @*/
4879371c9d4SSatish Balay PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed, void *ctx) {
48886d74e61SPeter Brune   PetscReal   angle = *(PetscReal *)linesearch->precheckctx;
48986d74e61SPeter Brune   Vec         Ylast;
49086d74e61SPeter Brune   PetscScalar dot;
49186d74e61SPeter Brune   PetscInt    iter;
49286d74e61SPeter Brune   PetscReal   ynorm, ylastnorm, theta, angle_radians;
49386d74e61SPeter Brune   SNES        snes;
49486d74e61SPeter Brune 
49586d74e61SPeter Brune   PetscFunctionBegin;
4969566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
4979566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)snes, "SNESLineSearchPreCheckPicard_Ylast", (PetscObject *)&Ylast));
49886d74e61SPeter Brune   if (!Ylast) {
4999566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(Y, &Ylast));
5009566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)snes, "SNESLineSearchPreCheckPicard_Ylast", (PetscObject)Ylast));
5019566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)Ylast));
50286d74e61SPeter Brune   }
5039566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &iter));
50486d74e61SPeter Brune   if (iter < 2) {
5059566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y, Ylast));
50686d74e61SPeter Brune     *changed = PETSC_FALSE;
50786d74e61SPeter Brune     PetscFunctionReturn(0);
50886d74e61SPeter Brune   }
50986d74e61SPeter Brune 
5109566063dSJacob Faibussowitsch   PetscCall(VecDot(Y, Ylast, &dot));
5119566063dSJacob Faibussowitsch   PetscCall(VecNorm(Y, NORM_2, &ynorm));
5129566063dSJacob Faibussowitsch   PetscCall(VecNorm(Ylast, NORM_2, &ylastnorm));
51386d74e61SPeter Brune   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
514255453a1SBarry Smith   theta         = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm), -1.0, 1.0));
51586d74e61SPeter Brune   angle_radians = angle * PETSC_PI / 180.;
51686d74e61SPeter Brune   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
51786d74e61SPeter Brune     /* Modify the step Y */
51886d74e61SPeter Brune     PetscReal alpha, ydiffnorm;
5199566063dSJacob Faibussowitsch     PetscCall(VecAXPY(Ylast, -1.0, Y));
5209566063dSJacob Faibussowitsch     PetscCall(VecNorm(Ylast, NORM_2, &ydiffnorm));
521f85e2ce2SBarry Smith     alpha = (ydiffnorm > .001 * ylastnorm) ? ylastnorm / ydiffnorm : 1000.0;
5229566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y, Ylast));
5239566063dSJacob Faibussowitsch     PetscCall(VecScale(Y, alpha));
5249566063dSJacob 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));
525a47ec85fSBarry Smith     *changed = PETSC_TRUE;
52686d74e61SPeter Brune   } else {
5279566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Angle %14.12e degrees exceeds threshold %14.12e, no correction applied\n", (double)(theta * 180. / PETSC_PI), (double)angle));
5289566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y, Ylast));
52986d74e61SPeter Brune     *changed = PETSC_FALSE;
530bf7f4e0aSPeter Brune   }
531bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
532bf7f4e0aSPeter Brune }
533bf7f4e0aSPeter Brune 
534f40b411bSPeter Brune /*@
535cd7522eaSPeter Brune    SNESLineSearchApply - Computes the line-search update.
536f40b411bSPeter Brune 
537*f6dfbefdSBarry Smith    Collective on linesearch
538f40b411bSPeter Brune 
539f40b411bSPeter Brune    Input Parameters:
5408e557f58SPeter Brune +  linesearch - The linesearch context
5418e557f58SPeter Brune -  Y - The search direction
542f40b411bSPeter Brune 
5436b867d5aSJose E. Roman    Input/Output Parameters:
5446b867d5aSJose E. Roman +  X - The current solution, on output the new solution
5456b867d5aSJose E. Roman .  F - The current function, on output the new function
5466b867d5aSJose E. Roman -  fnorm - The current norm, on output the new function norm
547f40b411bSPeter Brune 
548cd7522eaSPeter Brune    Options Database Keys:
5490b00b554SBarry Smith + -snes_linesearch_type - basic (or equivalently none), bt, l2, cp, nleqerr, shell
550dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches
5511fe24845SBarry Smith . -snes_linesearch_damping - The linesearch damping parameter, default is 1.0 (no damping)
5521fe24845SBarry Smith . -snes_linesearch_norms   - Turn on/off the linesearch norms computation (SNESLineSearchSetComputeNorms())
5533c7d6663SPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess
5543c7d6663SPeter Brune - -snes_linesearch_max_it - The number of iterations for iterative line searches
555cd7522eaSPeter Brune 
556cd7522eaSPeter Brune    Notes:
557*f6dfbefdSBarry Smith    This is typically called from within a `SNESSolve()` implementation in order to
558*f6dfbefdSBarry Smith    help with convergence of the nonlinear method.  Various `SNES` types use line searches
559cd7522eaSPeter Brune    in different ways, but the overarching theme is that a line search is used to determine
560cd7522eaSPeter Brune    an optimal damping parameter of a step at each iteration of the method.  Each
561*f6dfbefdSBarry Smith    application of the line search may invoke `SNESComputeFunction()` several times, and
562cd7522eaSPeter Brune    therefore may be fairly expensive.
563cd7522eaSPeter Brune 
564f40b411bSPeter Brune    Level: Intermediate
565f40b411bSPeter Brune 
566*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchPreCheck()`, `SNESLineSearchPostCheck()`, `SNESSolve()`, `SNESComputeFunction()`, `SNESLineSearchSetComputeNorms()`,
567db781477SPatrick Sanan           `SNESLineSearchType`, `SNESLineSearchSetType()`
568f40b411bSPeter Brune @*/
5699371c9d4SSatish Balay PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal *fnorm, Vec Y) {
570bf388a1fSBarry Smith   PetscFunctionBegin;
571f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
572bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
573bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
574064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(Y, VEC_CLASSID, 5);
575bf7f4e0aSPeter Brune 
576422a814eSBarry Smith   linesearch->result = SNES_LINESEARCH_SUCCEEDED;
577bf7f4e0aSPeter Brune 
578bf7f4e0aSPeter Brune   linesearch->vec_sol    = X;
579bf7f4e0aSPeter Brune   linesearch->vec_update = Y;
580bf7f4e0aSPeter Brune   linesearch->vec_func   = F;
581bf7f4e0aSPeter Brune 
5829566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchSetUp(linesearch));
583bf7f4e0aSPeter Brune 
584f5af7f23SKarl Rupp   if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */
585bf7f4e0aSPeter Brune 
586f5af7f23SKarl Rupp   if (fnorm) linesearch->fnorm = *fnorm;
58748a46eb9SPierre Jolivet   else PetscCall(VecNorm(F, NORM_2, &linesearch->fnorm));
588bf7f4e0aSPeter Brune 
5899566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(SNESLINESEARCH_Apply, linesearch, X, F, Y));
590bf7f4e0aSPeter Brune 
591dbbe0bcdSBarry Smith   PetscUseTypeMethod(linesearch, apply);
592bf7f4e0aSPeter Brune 
5939566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(SNESLINESEARCH_Apply, linesearch, X, F, Y));
594bf7f4e0aSPeter Brune 
595f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
596bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
597bf7f4e0aSPeter Brune }
598bf7f4e0aSPeter Brune 
599f40b411bSPeter Brune /*@
600f1c6b773SPeter Brune    SNESLineSearchDestroy - Destroys the line search instance.
601f40b411bSPeter Brune 
602*f6dfbefdSBarry Smith    Collective on linesearch
603f40b411bSPeter Brune 
604f40b411bSPeter Brune    Input Parameters:
6058e557f58SPeter Brune .  linesearch - The linesearch context
606f40b411bSPeter Brune 
60784238204SBarry Smith    Level: developer
608f40b411bSPeter Brune 
609*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchReset()`, `SNESDestroy()`
610f40b411bSPeter Brune @*/
6119371c9d4SSatish Balay PetscErrorCode SNESLineSearchDestroy(SNESLineSearch *linesearch) {
612bf7f4e0aSPeter Brune   PetscFunctionBegin;
613bf7f4e0aSPeter Brune   if (!*linesearch) PetscFunctionReturn(0);
614f1c6b773SPeter Brune   PetscValidHeaderSpecific((*linesearch), SNESLINESEARCH_CLASSID, 1);
6159371c9d4SSatish Balay   if (--((PetscObject)(*linesearch))->refct > 0) {
6169371c9d4SSatish Balay     *linesearch = NULL;
6179371c9d4SSatish Balay     PetscFunctionReturn(0);
6189371c9d4SSatish Balay   }
6199566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsViewOff((PetscObject)*linesearch));
6209566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchReset(*linesearch));
621f5af7f23SKarl Rupp   if ((*linesearch)->ops->destroy) (*linesearch)->ops->destroy(*linesearch);
6229566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&(*linesearch)->monitor));
6239566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchMonitorCancel((*linesearch)));
6249566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(linesearch));
625bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
626bf7f4e0aSPeter Brune }
627bf7f4e0aSPeter Brune 
628f40b411bSPeter Brune /*@
629dcf2fd19SBarry Smith    SNESLineSearchSetDefaultMonitor - Turns on/off printing useful information and debugging output about the line search.
630bf7f4e0aSPeter Brune 
631*f6dfbefdSBarry Smith    Logically Collective on linesearch
632*f6dfbefdSBarry Smith 
633bf7f4e0aSPeter Brune    Input Parameters:
634dcf2fd19SBarry Smith +  linesearch - the linesearch object
635*f6dfbefdSBarry Smith -  viewer - an `PETSCVIEWERASCII` `PetscViewer` or NULL to turn off monitor
636bf7f4e0aSPeter Brune 
637*f6dfbefdSBarry Smith    Options Database Key:
638dcf2fd19SBarry Smith .   -snes_linesearch_monitor [:filename] - enables the monitor
639bf7f4e0aSPeter Brune 
640bf7f4e0aSPeter Brune    Level: intermediate
641bf7f4e0aSPeter Brune 
642*f6dfbefdSBarry Smith    Developer Note:
643*f6dfbefdSBarry Smith    This monitor is implemented differently than the other line search monitors that are set with
644*f6dfbefdSBarry Smith    `SNESLineSearchMonitorSet()` since it is called in many locations of the line search routines to display aspects of the
645d12e167eSBarry Smith    line search that are not visible to the other monitors.
646dcf2fd19SBarry Smith 
647*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `PETSCVIEWERASCII`, `SNESGetLineSearch()`, `SNESLineSearchGetDefaultMonitor()`, `PetscViewer`, `SNESLineSearchSetMonitor()`,
648*f6dfbefdSBarry Smith           `SNESLineSearchMonitorSetFromOptions()`
649bf7f4e0aSPeter Brune @*/
6509371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetDefaultMonitor(SNESLineSearch linesearch, PetscViewer viewer) {
651bf7f4e0aSPeter Brune   PetscFunctionBegin;
6529566063dSJacob Faibussowitsch   if (viewer) PetscCall(PetscObjectReference((PetscObject)viewer));
6539566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&linesearch->monitor));
654dcf2fd19SBarry Smith   linesearch->monitor = viewer;
655bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
656bf7f4e0aSPeter Brune }
657bf7f4e0aSPeter Brune 
658f40b411bSPeter Brune /*@
659*f6dfbefdSBarry Smith    SNESLineSearchGetDefaultMonitor - Gets the `PetscViewer` instance for the line search monitor.
660*f6dfbefdSBarry Smith 
661*f6dfbefdSBarry Smith    Logically Collective on linsearch
6626a388c36SPeter Brune 
663f190f2fcSBarry Smith    Input Parameter:
6648e557f58SPeter Brune .  linesearch - linesearch context
665f40b411bSPeter Brune 
666f190f2fcSBarry Smith    Output Parameter:
6678e557f58SPeter Brune .  monitor - monitor context
668f40b411bSPeter Brune 
669f40b411bSPeter Brune    Level: intermediate
670f40b411bSPeter Brune 
671*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchSetDefaultMonitor()`, `PetscViewer`
672f40b411bSPeter Brune @*/
6739371c9d4SSatish Balay PetscErrorCode SNESLineSearchGetDefaultMonitor(SNESLineSearch linesearch, PetscViewer *monitor) {
6746a388c36SPeter Brune   PetscFunctionBegin;
675f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
6766a388c36SPeter Brune   *monitor = linesearch->monitor;
6776a388c36SPeter Brune   PetscFunctionReturn(0);
6786a388c36SPeter Brune }
6796a388c36SPeter Brune 
680dcf2fd19SBarry Smith /*@C
681dcf2fd19SBarry Smith    SNESLineSearchMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
682dcf2fd19SBarry Smith 
683*f6dfbefdSBarry Smith    Collective on ls
684dcf2fd19SBarry Smith 
685dcf2fd19SBarry Smith    Input Parameters:
686*f6dfbefdSBarry Smith +  ls - `SNESLineSearch` object you wish to monitor
687dcf2fd19SBarry Smith .  name - the monitor type one is seeking
688dcf2fd19SBarry Smith .  help - message indicating what monitoring is done
689dcf2fd19SBarry Smith .  manual - manual page for the monitor
690dcf2fd19SBarry Smith .  monitor - the monitor function
691*f6dfbefdSBarry 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`
692dcf2fd19SBarry Smith 
693dcf2fd19SBarry Smith    Level: developer
694dcf2fd19SBarry Smith 
695*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetMonitor()`, `PetscOptionsGetViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
696db781477SPatrick Sanan           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
697db781477SPatrick Sanan           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
698db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
699c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
700db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
701db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
702dcf2fd19SBarry Smith @*/
7039371c9d4SSatish Balay PetscErrorCode SNESLineSearchMonitorSetFromOptions(SNESLineSearch ls, const char name[], const char help[], const char manual[], PetscErrorCode (*monitor)(SNESLineSearch, PetscViewerAndFormat *), PetscErrorCode (*monitorsetup)(SNESLineSearch, PetscViewerAndFormat *)) {
704dcf2fd19SBarry Smith   PetscViewer       viewer;
705dcf2fd19SBarry Smith   PetscViewerFormat format;
706dcf2fd19SBarry Smith   PetscBool         flg;
707dcf2fd19SBarry Smith 
708dcf2fd19SBarry Smith   PetscFunctionBegin;
7099566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)ls), ((PetscObject)ls)->options, ((PetscObject)ls)->prefix, name, &viewer, &format, &flg));
710dcf2fd19SBarry Smith   if (flg) {
711d12e167eSBarry Smith     PetscViewerAndFormat *vf;
7129566063dSJacob Faibussowitsch     PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf));
7139566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)viewer));
7141baa6e33SBarry Smith     if (monitorsetup) PetscCall((*monitorsetup)(ls, vf));
7159566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchMonitorSet(ls, (PetscErrorCode(*)(SNESLineSearch, void *))monitor, vf, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy));
716dcf2fd19SBarry Smith   }
717dcf2fd19SBarry Smith   PetscFunctionReturn(0);
718dcf2fd19SBarry Smith }
719dcf2fd19SBarry Smith 
720f40b411bSPeter Brune /*@
721f1c6b773SPeter Brune    SNESLineSearchSetFromOptions - Sets options for the line search
722f40b411bSPeter Brune 
723*f6dfbefdSBarry Smith    Logically Collective on linesearch
724*f6dfbefdSBarry Smith 
725f899ff85SJose E. Roman    Input Parameter:
7268e557f58SPeter Brune .  linesearch - linesearch context
727f40b411bSPeter Brune 
728f40b411bSPeter Brune    Options Database Keys:
7290b00b554SBarry Smith + -snes_linesearch_type <type> - basic (or equivalently none), bt, l2, cp, nleqerr, shell
7303c7d6663SPeter Brune . -snes_linesearch_order <order> - 1, 2, 3.  Most types only support certain orders (bt supports 2 or 3)
731*f6dfbefdSBarry Smith . -snes_linesearch_norms   - Turn on/off the linesearch norms for the basic linesearch typem (`SNESLineSearchSetComputeNorms()`)
73271eef1aeSPeter Brune . -snes_linesearch_minlambda - The minimum step length
7331a9b3a06SPeter Brune . -snes_linesearch_maxstep - The maximum step size
7341a9b3a06SPeter Brune . -snes_linesearch_rtol - Relative tolerance for iterative line searches
7351a9b3a06SPeter Brune . -snes_linesearch_atol - Absolute tolerance for iterative line searches
7361a9b3a06SPeter Brune . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches
7371a9b3a06SPeter Brune . -snes_linesearch_max_it - The number of iterations for iterative line searches
738dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches
739dcf2fd19SBarry Smith . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine
7408e557f58SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter
741cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess.
7421a9b3a06SPeter Brune . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method
743d8d34be6SBarry Smith - -snes_linesearch_precheck_picard_angle - Angle used in Picard precheck method
744f40b411bSPeter Brune 
745f40b411bSPeter Brune    Level: intermediate
746f40b411bSPeter Brune 
747*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchCreate()`, `SNESLineSearchSetOrder()`, `SNESLineSearchSetType()`, `SNESLineSearchSetTolerances()`, `SNESLineSearchSetDamping()`, `SNESLineSearchPreCheckPicard()`,
748db781477SPatrick Sanan           `SNESLineSearchType`, `SNESLineSearchSetComputeNorms()`
749f40b411bSPeter Brune @*/
7509371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch) {
7511a4f838cSPeter Brune   const char *deft = SNESLINESEARCHBASIC;
752bf7f4e0aSPeter Brune   char        type[256];
753bf7f4e0aSPeter Brune   PetscBool   flg, set;
754dcf2fd19SBarry Smith   PetscViewer viewer;
755bf388a1fSBarry Smith 
756bf7f4e0aSPeter Brune   PetscFunctionBegin;
7579566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchRegisterAll());
758bf7f4e0aSPeter Brune 
759d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)linesearch);
760f5af7f23SKarl Rupp   if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name;
7619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-snes_linesearch_type", "Linesearch type", "SNESLineSearchSetType", SNESLineSearchList, deft, type, 256, &flg));
762bf7f4e0aSPeter Brune   if (flg) {
7639566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetType(linesearch, type));
764bf7f4e0aSPeter Brune   } else if (!((PetscObject)linesearch)->type_name) {
7659566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetType(linesearch, deft));
766bf7f4e0aSPeter Brune   }
767bf7f4e0aSPeter Brune 
7689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)linesearch), ((PetscObject)linesearch)->options, ((PetscObject)linesearch)->prefix, "-snes_linesearch_monitor", &viewer, NULL, &set));
769dcf2fd19SBarry Smith   if (set) {
7709566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetDefaultMonitor(linesearch, viewer));
7719566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
772dcf2fd19SBarry Smith   }
7739566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchMonitorSetFromOptions(linesearch, "-snes_linesearch_monitor_solution_update", "View correction at each iteration", "SNESLineSearchMonitorSolutionUpdate", SNESLineSearchMonitorSolutionUpdate, NULL));
774bf7f4e0aSPeter Brune 
7751a9b3a06SPeter Brune   /* tolerances */
7769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_minlambda", "Minimum step length", "SNESLineSearchSetTolerances", linesearch->steptol, &linesearch->steptol, NULL));
7779566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_maxstep", "Maximum step size", "SNESLineSearchSetTolerances", linesearch->maxstep, &linesearch->maxstep, NULL));
7789566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_rtol", "Relative tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->rtol, &linesearch->rtol, NULL));
7799566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_atol", "Absolute tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->atol, &linesearch->atol, NULL));
7809566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_ltol", "Change in lambda tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->ltol, &linesearch->ltol, NULL));
7819566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_linesearch_max_it", "Maximum iterations for iterative line searches", "SNESLineSearchSetTolerances", linesearch->max_its, &linesearch->max_its, NULL));
7827a35526eSPeter Brune 
7831a9b3a06SPeter Brune   /* damping parameters */
7849566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_damping", "Line search damping and initial step guess", "SNESLineSearchSetDamping", linesearch->damping, &linesearch->damping, NULL));
7851a9b3a06SPeter Brune 
7869566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_linesearch_keeplambda", "Use previous lambda as damping", "SNESLineSearchSetKeepLambda", linesearch->keeplambda, &linesearch->keeplambda, NULL));
7871a9b3a06SPeter Brune 
7881a9b3a06SPeter Brune   /* precheck */
7899566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_linesearch_precheck_picard", "Use a correction that sometimes improves convergence of Picard iteration", "SNESLineSearchPreCheckPicard", flg, &flg, &set));
7907a35526eSPeter Brune   if (set) {
7917a35526eSPeter Brune     if (flg) {
7927a35526eSPeter Brune       linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */
793f5af7f23SKarl Rupp 
794d0609cedSBarry 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));
7959566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchSetPreCheck(linesearch, SNESLineSearchPreCheckPicard, &linesearch->precheck_picard_angle));
7967a35526eSPeter Brune     } else {
7979566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchSetPreCheck(linesearch, NULL, NULL));
7987a35526eSPeter Brune     }
7997a35526eSPeter Brune   }
8009566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_linesearch_order", "Order of approximation used in the line search", "SNESLineSearchSetOrder", linesearch->order, &linesearch->order, NULL));
8019566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_linesearch_norms", "Compute final norms in line search", "SNESLineSearchSetComputeNorms", linesearch->norms, &linesearch->norms, NULL));
8027a35526eSPeter Brune 
803dbbe0bcdSBarry Smith   PetscTryTypeMethod(linesearch, setfromoptions, PetscOptionsObject);
8045a9b6599SPeter Brune 
805dbbe0bcdSBarry Smith   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)linesearch, PetscOptionsObject));
806d0609cedSBarry Smith   PetscOptionsEnd();
807bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
808bf7f4e0aSPeter Brune }
809bf7f4e0aSPeter Brune 
810f40b411bSPeter Brune /*@
811f190f2fcSBarry Smith    SNESLineSearchView - Prints useful information about the line search
812f40b411bSPeter Brune 
813f40b411bSPeter Brune    Input Parameters:
8148e557f58SPeter Brune .  linesearch - linesearch context
815f40b411bSPeter Brune 
816*f6dfbefdSBarry Smith    Logically Collective on linesearch
817f40b411bSPeter Brune 
818f40b411bSPeter Brune    Level: intermediate
819f40b411bSPeter Brune 
820*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `PetscViewer`, `SNESLineSearchCreate()`
821f40b411bSPeter Brune @*/
8229371c9d4SSatish Balay PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer) {
8237f1410a3SPeter Brune   PetscBool iascii;
824bf388a1fSBarry Smith 
825bf7f4e0aSPeter Brune   PetscFunctionBegin;
8267f1410a3SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
82748a46eb9SPierre Jolivet   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch), &viewer));
8287f1410a3SPeter Brune   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
8297f1410a3SPeter Brune   PetscCheckSameComm(linesearch, 1, viewer, 2);
830f40b411bSPeter Brune 
8319566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
8327f1410a3SPeter Brune   if (iascii) {
8339566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)linesearch, viewer));
8349566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
835dbbe0bcdSBarry Smith     PetscTryTypeMethod(linesearch, view, viewer);
8369566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
8379566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep, (double)linesearch->steptol));
8389566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol, (double)linesearch->atol, (double)linesearch->ltol));
83963a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  maximum iterations=%" PetscInt_FMT "\n", linesearch->max_its));
8406b2b7091SBarry Smith     if (linesearch->ops->precheck) {
8416b2b7091SBarry Smith       if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) {
84263a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  using precheck step to speed up Picard convergence\n"));
8437f1410a3SPeter Brune       } else {
84463a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  using user-defined precheck step\n"));
8457f1410a3SPeter Brune       }
8467f1410a3SPeter Brune     }
84748a46eb9SPierre Jolivet     if (linesearch->ops->postcheck) PetscCall(PetscViewerASCIIPrintf(viewer, "  using user-defined postcheck step\n"));
8487f1410a3SPeter Brune   }
849bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
850bf7f4e0aSPeter Brune }
851bf7f4e0aSPeter Brune 
852ea5d4fccSPeter Brune /*@C
853a80ff896SJed Brown    SNESLineSearchGetType - Gets the linesearch type
854a80ff896SJed Brown 
855*f6dfbefdSBarry Smith    Logically Collective on linesearch
856a80ff896SJed Brown 
857a80ff896SJed Brown    Input Parameters:
858a80ff896SJed Brown .  linesearch - linesearch context
859a80ff896SJed Brown 
860a80ff896SJed Brown    Output Parameters:
861a80ff896SJed Brown -  type - The type of line search, or NULL if not set
862a80ff896SJed Brown 
863a80ff896SJed Brown    Level: intermediate
864a80ff896SJed Brown 
865*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchType`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchSetType()`
866a80ff896SJed Brown @*/
8679371c9d4SSatish Balay PetscErrorCode SNESLineSearchGetType(SNESLineSearch linesearch, SNESLineSearchType *type) {
868a80ff896SJed Brown   PetscFunctionBegin;
869a80ff896SJed Brown   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
870dadcf809SJacob Faibussowitsch   PetscValidPointer(type, 2);
871a80ff896SJed Brown   *type = ((PetscObject)linesearch)->type_name;
872a80ff896SJed Brown   PetscFunctionReturn(0);
873a80ff896SJed Brown }
874a80ff896SJed Brown 
875a80ff896SJed Brown /*@C
876f1c6b773SPeter Brune    SNESLineSearchSetType - Sets the linesearch type
877f40b411bSPeter Brune 
878*f6dfbefdSBarry Smith    Logically Collective on linesearch
879f190f2fcSBarry Smith 
880f40b411bSPeter Brune    Input Parameters:
8818e557f58SPeter Brune +  linesearch - linesearch context
882f40b411bSPeter Brune -  type - The type of line search to be used
883f40b411bSPeter Brune 
884cd7522eaSPeter Brune    Available Types:
885*f6dfbefdSBarry Smith +  `SNESLINESEARCHBASIC` - (or equivalently `SNESLINESEARCHNONE`) Simple damping line search, defaults to using the full Newton step
886*f6dfbefdSBarry Smith .  `SNESLINESEARCHBT` - Backtracking line search over the L2 norm of the function
887*f6dfbefdSBarry Smith .  `SNESLINESEARCHL2` - Secant line search over the L2 norm of the function
888*f6dfbefdSBarry Smith .  `SNESLINESEARCHCP` - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x)
889*f6dfbefdSBarry Smith .  `SNESLINESEARCHNLEQERR` - Affine-covariant error-oriented linesearch
890*f6dfbefdSBarry Smith -  `SNESLINESEARCHSHELL` - User provided `SNESLineSearch` implementation
8911fe24845SBarry Smith 
8921fe24845SBarry Smith    Options Database:
8930b00b554SBarry Smith .  -snes_linesearch_type <type> - basic (or equivalently none), bt, l2, cp, nleqerr, shell
894cd7522eaSPeter Brune 
895f40b411bSPeter Brune    Level: intermediate
896f40b411bSPeter Brune 
897*f6dfbefdSBarry Smith .seealso:  `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchType`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchGetType()`
898f40b411bSPeter Brune @*/
8999371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type) {
900bf7f4e0aSPeter Brune   PetscBool match;
9015f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(SNESLineSearch);
902bf7f4e0aSPeter Brune 
903bf7f4e0aSPeter Brune   PetscFunctionBegin;
904f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
905bf7f4e0aSPeter Brune   PetscValidCharPointer(type, 2);
906bf7f4e0aSPeter Brune 
9079566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)linesearch, type, &match));
908bf7f4e0aSPeter Brune   if (match) PetscFunctionReturn(0);
909bf7f4e0aSPeter Brune 
9109566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(SNESLineSearchList, type, &r));
9115f80ce2aSJacob Faibussowitsch   PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested Line Search type %s", type);
912bf7f4e0aSPeter Brune   /* Destroy the previous private linesearch context */
913bf7f4e0aSPeter Brune   if (linesearch->ops->destroy) {
9149566063dSJacob Faibussowitsch     PetscCall((*(linesearch)->ops->destroy)(linesearch));
9150298fd71SBarry Smith     linesearch->ops->destroy = NULL;
916bf7f4e0aSPeter Brune   }
917f1c6b773SPeter Brune   /* Reinitialize function pointers in SNESLineSearchOps structure */
9189e5d0892SLisandro Dalcin   linesearch->ops->apply          = NULL;
9199e5d0892SLisandro Dalcin   linesearch->ops->view           = NULL;
9209e5d0892SLisandro Dalcin   linesearch->ops->setfromoptions = NULL;
9219e5d0892SLisandro Dalcin   linesearch->ops->destroy        = NULL;
922bf7f4e0aSPeter Brune 
9239566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)linesearch, type));
9249566063dSJacob Faibussowitsch   PetscCall((*r)(linesearch));
925bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
926bf7f4e0aSPeter Brune }
927bf7f4e0aSPeter Brune 
928f40b411bSPeter Brune /*@
929*f6dfbefdSBarry Smith    SNESLineSearchSetSNES - Sets the `SNES` for the linesearch for function evaluation.
930f40b411bSPeter Brune 
931f40b411bSPeter Brune    Input Parameters:
9328e557f58SPeter Brune +  linesearch - linesearch context
933f40b411bSPeter Brune -  snes - The snes instance
934f40b411bSPeter Brune 
93578bcb3b5SPeter Brune    Level: developer
93678bcb3b5SPeter Brune 
937*f6dfbefdSBarry Smith    Note:
938f190f2fcSBarry Smith    This happens automatically when the line search is obtained/created with
939*f6dfbefdSBarry Smith    `SNESGetLineSearch()`.  This routine is therefore mainly called within `SNES`
94078bcb3b5SPeter Brune    implementations.
941f40b411bSPeter Brune 
9428141a3b9SPeter Brune    Level: developer
943f40b411bSPeter Brune 
944*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetSNES()`, `SNESLineSearchSetVecs()`, `SNES`
945f40b411bSPeter Brune @*/
9469371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes) {
947bf7f4e0aSPeter Brune   PetscFunctionBegin;
948f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
949bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 2);
950bf7f4e0aSPeter Brune   linesearch->snes = snes;
951bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
952bf7f4e0aSPeter Brune }
953bf7f4e0aSPeter Brune 
954f40b411bSPeter Brune /*@
955*f6dfbefdSBarry Smith    SNESLineSearchGetSNES - Gets the `SNES` instance associated with the line search.
956*f6dfbefdSBarry Smith    Having an associated `SNES` is necessary because most line search implementations must be able to
957*f6dfbefdSBarry Smith    evaluate the function using `SNESComputeFunction()` for the associated `SNES`.  This routine
958*f6dfbefdSBarry Smith    is used in the line search implementations when one must get this associated `SNES` instance.
959*f6dfbefdSBarry Smith 
960*f6dfbefdSBarry Smith    Not Collective
961f40b411bSPeter Brune 
962f40b411bSPeter Brune    Input Parameters:
9638e557f58SPeter Brune .  linesearch - linesearch context
964f40b411bSPeter Brune 
965f40b411bSPeter Brune    Output Parameters:
966f40b411bSPeter Brune .  snes - The snes instance
967f40b411bSPeter Brune 
9688141a3b9SPeter Brune    Level: developer
969f40b411bSPeter Brune 
970*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESType`, `SNESLineSearchGetSNES()`, `SNESLineSearchSetVecs()`, `SNES`
971f40b411bSPeter Brune @*/
9729371c9d4SSatish Balay PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes) {
973bf7f4e0aSPeter Brune   PetscFunctionBegin;
974f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
9756a388c36SPeter Brune   PetscValidPointer(snes, 2);
976bf7f4e0aSPeter Brune   *snes = linesearch->snes;
977bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
978bf7f4e0aSPeter Brune }
979bf7f4e0aSPeter Brune 
980f40b411bSPeter Brune /*@
981f1c6b773SPeter Brune    SNESLineSearchGetLambda - Gets the last linesearch steplength discovered.
982f40b411bSPeter Brune 
983*f6dfbefdSBarry Smith    Not Collective
984*f6dfbefdSBarry Smith 
985f40b411bSPeter Brune    Input Parameters:
9868e557f58SPeter Brune .  linesearch - linesearch context
987f40b411bSPeter Brune 
988f40b411bSPeter Brune    Output Parameters:
989*f6dfbefdSBarry Smith .  lambda - The last steplength computed during `SNESLineSearchApply()`
990f40b411bSPeter Brune 
99178bcb3b5SPeter Brune    Level: advanced
99278bcb3b5SPeter Brune 
993*f6dfbefdSBarry Smith    Note:
9948e557f58SPeter Brune    This is useful in methods where the solver is ill-scaled and
99578bcb3b5SPeter Brune    requires some adaptive notion of the difference in scale between the
996*f6dfbefdSBarry Smith    solution and the function.  For instance, `SNESQN` may be scaled by the
99778bcb3b5SPeter Brune    line search lambda using the argument -snes_qn_scaling ls.
99878bcb3b5SPeter Brune 
999*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetLambda()`, `SNESLineSearchGetDamping()`, `SNESLineSearchApply()`
1000f40b411bSPeter Brune @*/
10019371c9d4SSatish Balay PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch linesearch, PetscReal *lambda) {
10026a388c36SPeter Brune   PetscFunctionBegin;
1003f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1004534a8f05SLisandro Dalcin   PetscValidRealPointer(lambda, 2);
10056a388c36SPeter Brune   *lambda = linesearch->lambda;
10066a388c36SPeter Brune   PetscFunctionReturn(0);
10076a388c36SPeter Brune }
10086a388c36SPeter Brune 
1009f40b411bSPeter Brune /*@
1010*f6dfbefdSBarry Smith    SNESLineSearchSetLambda - Sets the linesearch steplength
1011f40b411bSPeter Brune 
1012f40b411bSPeter Brune    Input Parameters:
10138e557f58SPeter Brune +  linesearch - linesearch context
1014f40b411bSPeter Brune -  lambda - The last steplength.
1015f40b411bSPeter Brune 
1016*f6dfbefdSBarry Smith    Note:
1017*f6dfbefdSBarry Smith    This routine is typically used within implementations of `SNESLineSearchApply()`
1018*f6dfbefdSBarry Smith    to set the final steplength.  This routine (and `SNESLineSearchGetLambda()`) were
1019cd7522eaSPeter Brune    added in order to facilitate Quasi-Newton methods that use the previous steplength
1020cd7522eaSPeter Brune    as an inner scaling parameter.
1021cd7522eaSPeter Brune 
102278bcb3b5SPeter Brune    Level: advanced
1023f40b411bSPeter Brune 
1024*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetLambda()`
1025f40b411bSPeter Brune @*/
10269371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda) {
10276a388c36SPeter Brune   PetscFunctionBegin;
1028f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
10296a388c36SPeter Brune   linesearch->lambda = lambda;
10306a388c36SPeter Brune   PetscFunctionReturn(0);
10316a388c36SPeter Brune }
10326a388c36SPeter Brune 
1033f40b411bSPeter Brune /*@
10343c7d6663SPeter Brune    SNESLineSearchGetTolerances - Gets the tolerances for the linesearch.  These include
103578bcb3b5SPeter Brune    tolerances for the relative and absolute change in the function norm, the change
103678bcb3b5SPeter Brune    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
103778bcb3b5SPeter Brune    and the maximum number of iterations the line search procedure may take.
1038f40b411bSPeter Brune 
1039*f6dfbefdSBarry Smith    Not Collective
1040*f6dfbefdSBarry Smith 
1041f899ff85SJose E. Roman    Input Parameter:
10428e557f58SPeter Brune .  linesearch - linesearch context
1043f40b411bSPeter Brune 
1044f40b411bSPeter Brune    Output Parameters:
1045516fe3c3SPeter Brune +  steptol - The minimum steplength
10466cc8e53bSPeter Brune .  maxstep - The maximum steplength
1047516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
1048516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
1049516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
1050516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
1051f40b411bSPeter Brune 
105278bcb3b5SPeter Brune    Level: intermediate
105378bcb3b5SPeter Brune 
1054*f6dfbefdSBarry Smith    Note:
105578bcb3b5SPeter Brune    Different line searches may implement these parameters slightly differently as
10563c7d6663SPeter Brune    the type requires.
1057516fe3c3SPeter Brune 
1058*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetTolerances()`
1059f40b411bSPeter Brune @*/
10609371c9d4SSatish Balay PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch linesearch, PetscReal *steptol, PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its) {
10616a388c36SPeter Brune   PetscFunctionBegin;
1062f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1063516fe3c3SPeter Brune   if (steptol) {
1064534a8f05SLisandro Dalcin     PetscValidRealPointer(steptol, 2);
10656a388c36SPeter Brune     *steptol = linesearch->steptol;
1066516fe3c3SPeter Brune   }
1067516fe3c3SPeter Brune   if (maxstep) {
1068534a8f05SLisandro Dalcin     PetscValidRealPointer(maxstep, 3);
1069516fe3c3SPeter Brune     *maxstep = linesearch->maxstep;
1070516fe3c3SPeter Brune   }
1071516fe3c3SPeter Brune   if (rtol) {
1072534a8f05SLisandro Dalcin     PetscValidRealPointer(rtol, 4);
1073516fe3c3SPeter Brune     *rtol = linesearch->rtol;
1074516fe3c3SPeter Brune   }
1075516fe3c3SPeter Brune   if (atol) {
1076534a8f05SLisandro Dalcin     PetscValidRealPointer(atol, 5);
1077516fe3c3SPeter Brune     *atol = linesearch->atol;
1078516fe3c3SPeter Brune   }
1079516fe3c3SPeter Brune   if (ltol) {
1080534a8f05SLisandro Dalcin     PetscValidRealPointer(ltol, 6);
1081516fe3c3SPeter Brune     *ltol = linesearch->ltol;
1082516fe3c3SPeter Brune   }
1083516fe3c3SPeter Brune   if (max_its) {
1084534a8f05SLisandro Dalcin     PetscValidIntPointer(max_its, 7);
1085516fe3c3SPeter Brune     *max_its = linesearch->max_its;
1086516fe3c3SPeter Brune   }
10876a388c36SPeter Brune   PetscFunctionReturn(0);
10886a388c36SPeter Brune }
10896a388c36SPeter Brune 
1090f40b411bSPeter Brune /*@
10913c7d6663SPeter Brune    SNESLineSearchSetTolerances -  Gets the tolerances for the linesearch.  These include
109278bcb3b5SPeter Brune    tolerances for the relative and absolute change in the function norm, the change
109378bcb3b5SPeter Brune    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
109478bcb3b5SPeter Brune    and the maximum number of iterations the line search procedure may take.
1095f40b411bSPeter Brune 
1096*f6dfbefdSBarry Smith    Collective
1097*f6dfbefdSBarry Smith 
1098f40b411bSPeter Brune    Input Parameters:
10998e557f58SPeter Brune +  linesearch - linesearch context
1100516fe3c3SPeter Brune .  steptol - The minimum steplength
11016cc8e53bSPeter Brune .  maxstep - The maximum steplength
1102516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
1103516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
1104516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
1105516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
1106f40b411bSPeter Brune 
1107*f6dfbefdSBarry Smith    Note:
1108*f6dfbefdSBarry Smith    The user may choose to not set any of the tolerances using `PETSC_DEFAULT` in
110978bcb3b5SPeter Brune    place of an argument.
1110f40b411bSPeter Brune 
111178bcb3b5SPeter Brune    Level: intermediate
1112516fe3c3SPeter Brune 
1113*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetTolerances()`
1114f40b411bSPeter Brune @*/
11159371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch linesearch, PetscReal steptol, PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its) {
11166a388c36SPeter Brune   PetscFunctionBegin;
1117f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1118d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, steptol, 2);
1119d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, maxstep, 3);
1120d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, rtol, 4);
1121d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, atol, 5);
1122d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, ltol, 6);
1123d3952184SSatish Balay   PetscValidLogicalCollectiveInt(linesearch, max_its, 7);
1124d3952184SSatish Balay 
1125d3952184SSatish Balay   if (steptol != PETSC_DEFAULT) {
11265f80ce2aSJacob Faibussowitsch     PetscCheck(steptol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Minimum step length %14.12e must be non-negative", (double)steptol);
11276a388c36SPeter Brune     linesearch->steptol = steptol;
1128d3952184SSatish Balay   }
1129d3952184SSatish Balay 
1130d3952184SSatish Balay   if (maxstep != PETSC_DEFAULT) {
11315f80ce2aSJacob Faibussowitsch     PetscCheck(maxstep >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum step length %14.12e must be non-negative", (double)maxstep);
1132516fe3c3SPeter Brune     linesearch->maxstep = maxstep;
1133d3952184SSatish Balay   }
1134d3952184SSatish Balay 
1135d3952184SSatish Balay   if (rtol != PETSC_DEFAULT) {
11362061ca32SJunchao 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);
1137516fe3c3SPeter Brune     linesearch->rtol = rtol;
1138d3952184SSatish Balay   }
1139d3952184SSatish Balay 
1140d3952184SSatish Balay   if (atol != PETSC_DEFAULT) {
11415f80ce2aSJacob Faibussowitsch     PetscCheck(atol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Absolute tolerance %14.12e must be non-negative", (double)atol);
1142516fe3c3SPeter Brune     linesearch->atol = atol;
1143d3952184SSatish Balay   }
1144d3952184SSatish Balay 
1145d3952184SSatish Balay   if (ltol != PETSC_DEFAULT) {
11465f80ce2aSJacob Faibussowitsch     PetscCheck(ltol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Lambda tolerance %14.12e must be non-negative", (double)ltol);
1147516fe3c3SPeter Brune     linesearch->ltol = ltol;
1148d3952184SSatish Balay   }
1149d3952184SSatish Balay 
1150d3952184SSatish Balay   if (max_its != PETSC_DEFAULT) {
115163a3b9bcSJacob Faibussowitsch     PetscCheck(max_its >= 0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of iterations %" PetscInt_FMT " must be non-negative", max_its);
1152516fe3c3SPeter Brune     linesearch->max_its = max_its;
1153d3952184SSatish Balay   }
11546a388c36SPeter Brune   PetscFunctionReturn(0);
11556a388c36SPeter Brune }
11566a388c36SPeter Brune 
1157f40b411bSPeter Brune /*@
1158f1c6b773SPeter Brune    SNESLineSearchGetDamping - Gets the line search damping parameter.
1159f40b411bSPeter Brune 
1160f40b411bSPeter Brune    Input Parameters:
11618e557f58SPeter Brune .  linesearch - linesearch context
1162f40b411bSPeter Brune 
1163f40b411bSPeter Brune    Output Parameters:
11648e557f58SPeter Brune .  damping - The damping parameter
1165f40b411bSPeter Brune 
116678bcb3b5SPeter Brune    Level: advanced
1167f40b411bSPeter Brune 
1168db781477SPatrick Sanan .seealso: `SNESLineSearchGetStepTolerance()`, `SNESQN`
1169f40b411bSPeter Brune @*/
1170f40b411bSPeter Brune 
11719371c9d4SSatish Balay PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch linesearch, PetscReal *damping) {
11726a388c36SPeter Brune   PetscFunctionBegin;
1173f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1174534a8f05SLisandro Dalcin   PetscValidRealPointer(damping, 2);
11756a388c36SPeter Brune   *damping = linesearch->damping;
11766a388c36SPeter Brune   PetscFunctionReturn(0);
11776a388c36SPeter Brune }
11786a388c36SPeter Brune 
1179f40b411bSPeter Brune /*@
1180fd292e60Sprj-    SNESLineSearchSetDamping - Sets the line search damping parameter.
1181f40b411bSPeter Brune 
1182f40b411bSPeter Brune    Input Parameters:
118303fd4120SBarry Smith +  linesearch - linesearch context
118403fd4120SBarry Smith -  damping - The damping parameter
1185f40b411bSPeter Brune 
1186*f6dfbefdSBarry Smith    Options Database Key:
1187*f6dfbefdSBarry Smith .   -snes_linesearch_damping <damping> - the damping value
1188*f6dfbefdSBarry Smith 
1189f40b411bSPeter Brune    Level: intermediate
1190f40b411bSPeter Brune 
1191*f6dfbefdSBarry Smith    Note:
1192*f6dfbefdSBarry Smith    The `SNESLINESEARCHNONE` line search merely takes the update step scaled by the damping parameter.
1193cd7522eaSPeter Brune    The use of the damping parameter in the l2 and cp line searches is much more subtle;
119478bcb3b5SPeter Brune    it is used as a starting point in calculating the secant step. However, the eventual
1195cd7522eaSPeter Brune    step may be of greater length than the damping parameter.  In the bt line search it is
1196cd7522eaSPeter Brune    used as the maximum possible step length, as the bt line search only backtracks.
1197cd7522eaSPeter Brune 
1198*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetDamping()`
1199f40b411bSPeter Brune @*/
12009371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch linesearch, PetscReal damping) {
12016a388c36SPeter Brune   PetscFunctionBegin;
1202f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
12036a388c36SPeter Brune   linesearch->damping = damping;
12046a388c36SPeter Brune   PetscFunctionReturn(0);
12056a388c36SPeter Brune }
12066a388c36SPeter Brune 
120759405d5eSPeter Brune /*@
120859405d5eSPeter Brune    SNESLineSearchGetOrder - Gets the line search approximation order.
120959405d5eSPeter Brune 
1210*f6dfbefdSBarry Smith    Input Parameter:
121178bcb3b5SPeter Brune .  linesearch - linesearch context
121259405d5eSPeter Brune 
1213*f6dfbefdSBarry Smith    Output Parameter:
12148e557f58SPeter Brune .  order - The order
121559405d5eSPeter Brune 
121678bcb3b5SPeter Brune    Possible Values for order:
1217*f6dfbefdSBarry Smith +  1 or `SNES_LINESEARCH_ORDER_LINEAR` - linear order
1218*f6dfbefdSBarry Smith .  2 or `SNES_LINESEARCH_ORDER_QUADRATIC` - quadratic order
1219*f6dfbefdSBarry Smith -  3 or `SNES_LINESEARCH_ORDER_CUBIC` - cubic order
122078bcb3b5SPeter Brune 
122159405d5eSPeter Brune    Level: intermediate
122259405d5eSPeter Brune 
1223*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetOrder()`
122459405d5eSPeter Brune @*/
122559405d5eSPeter Brune 
12269371c9d4SSatish Balay PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch linesearch, PetscInt *order) {
122759405d5eSPeter Brune   PetscFunctionBegin;
122859405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1229534a8f05SLisandro Dalcin   PetscValidIntPointer(order, 2);
123059405d5eSPeter Brune   *order = linesearch->order;
123159405d5eSPeter Brune   PetscFunctionReturn(0);
123259405d5eSPeter Brune }
123359405d5eSPeter Brune 
123459405d5eSPeter Brune /*@
12351f8196a2SJed Brown    SNESLineSearchSetOrder - Sets the maximum order of the polynomial fit used in the line search
123659405d5eSPeter Brune 
123759405d5eSPeter Brune    Input Parameters:
1238a2b725a8SWilliam Gropp +  linesearch - linesearch context
1239a2b725a8SWilliam Gropp -  order - The damping parameter
124059405d5eSPeter Brune 
124159405d5eSPeter Brune    Level: intermediate
124259405d5eSPeter Brune 
124378bcb3b5SPeter Brune    Possible Values for order:
1244*f6dfbefdSBarry Smith +  1 or `SNES_LINESEARCH_ORDER_LINEAR` - linear order
1245*f6dfbefdSBarry Smith .  2 or `SNES_LINESEARCH_ORDER_QUADRATIC` - quadratic order
1246*f6dfbefdSBarry Smith -  3 or `SNES_LINESEARCH_ORDER_CUBIC` - cubic order
124778bcb3b5SPeter Brune 
124859405d5eSPeter Brune    Notes:
124959405d5eSPeter Brune    Variable orders are supported by the following line searches:
125078bcb3b5SPeter Brune +  bt - cubic and quadratic
125178bcb3b5SPeter Brune -  cp - linear and quadratic
125259405d5eSPeter Brune 
1253*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetOrder()`, `SNESLineSearchSetDamping()`
125459405d5eSPeter Brune @*/
12559371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch linesearch, PetscInt order) {
125659405d5eSPeter Brune   PetscFunctionBegin;
125759405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
125859405d5eSPeter Brune   linesearch->order = order;
125959405d5eSPeter Brune   PetscFunctionReturn(0);
126059405d5eSPeter Brune }
126159405d5eSPeter Brune 
1262f40b411bSPeter Brune /*@
1263f1c6b773SPeter Brune    SNESLineSearchGetNorms - Gets the norms for for X, Y, and F.
1264f40b411bSPeter Brune 
1265*f6dfbefdSBarry Smith    Not Collective
1266*f6dfbefdSBarry Smith 
1267f899ff85SJose E. Roman    Input Parameter:
126878bcb3b5SPeter Brune .  linesearch - linesearch context
1269f40b411bSPeter Brune 
1270f40b411bSPeter Brune    Output Parameters:
1271f40b411bSPeter Brune +  xnorm - The norm of the current solution
1272f40b411bSPeter Brune .  fnorm - The norm of the current function
1273f40b411bSPeter Brune -  ynorm - The norm of the current update
1274f40b411bSPeter Brune 
127578bcb3b5SPeter Brune    Level: developer
1276f40b411bSPeter Brune 
1277*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetNorms()` `SNESLineSearchGetVecs()`
1278f40b411bSPeter Brune @*/
12799371c9d4SSatish Balay PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal *xnorm, PetscReal *fnorm, PetscReal *ynorm) {
1280bf7f4e0aSPeter Brune   PetscFunctionBegin;
1281f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1282f5af7f23SKarl Rupp   if (xnorm) *xnorm = linesearch->xnorm;
1283f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
1284f5af7f23SKarl Rupp   if (ynorm) *ynorm = linesearch->ynorm;
1285bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1286bf7f4e0aSPeter Brune }
1287bf7f4e0aSPeter Brune 
1288f40b411bSPeter Brune /*@
1289f1c6b773SPeter Brune    SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F.
1290f40b411bSPeter Brune 
1291*f6dfbefdSBarry Smith    Collective on linesearch
1292*f6dfbefdSBarry Smith 
1293f40b411bSPeter Brune    Input Parameters:
129478bcb3b5SPeter Brune +  linesearch - linesearch context
1295f40b411bSPeter Brune .  xnorm - The norm of the current solution
1296f40b411bSPeter Brune .  fnorm - The norm of the current function
1297f40b411bSPeter Brune -  ynorm - The norm of the current update
1298f40b411bSPeter Brune 
1299*f6dfbefdSBarry Smith    Level: developer
1300f40b411bSPeter Brune 
1301*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()`
1302f40b411bSPeter Brune @*/
13039371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm) {
13046a388c36SPeter Brune   PetscFunctionBegin;
1305f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
13066a388c36SPeter Brune   linesearch->xnorm = xnorm;
13076a388c36SPeter Brune   linesearch->fnorm = fnorm;
13086a388c36SPeter Brune   linesearch->ynorm = ynorm;
13096a388c36SPeter Brune   PetscFunctionReturn(0);
13106a388c36SPeter Brune }
13116a388c36SPeter Brune 
1312f40b411bSPeter Brune /*@
1313f1c6b773SPeter Brune    SNESLineSearchComputeNorms - Computes the norms of X, F, and Y.
1314f40b411bSPeter Brune 
1315f40b411bSPeter Brune    Input Parameters:
131678bcb3b5SPeter Brune .  linesearch - linesearch context
1317f40b411bSPeter Brune 
1318f40b411bSPeter Brune    Options Database Keys:
13198e557f58SPeter Brune .   -snes_linesearch_norms - turn norm computation on or off
1320f40b411bSPeter Brune 
1321f40b411bSPeter Brune    Level: intermediate
1322f40b411bSPeter Brune 
1323*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetNorms`, `SNESLineSearchSetNorms()`, `SNESLineSearchSetComputeNorms()`
1324f40b411bSPeter Brune @*/
13259371c9d4SSatish Balay PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch) {
13269bd66eb0SPeter Brune   SNES snes;
1327bf388a1fSBarry Smith 
13286a388c36SPeter Brune   PetscFunctionBegin;
13296a388c36SPeter Brune   if (linesearch->norms) {
13309bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
13319566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
13329566063dSJacob Faibussowitsch       PetscCall(VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
13339566063dSJacob Faibussowitsch       PetscCall(VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm));
13349566063dSJacob Faibussowitsch       PetscCall((*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm));
13359bd66eb0SPeter Brune     } else {
13369566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(linesearch->vec_func, NORM_2, &linesearch->fnorm));
13379566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
13389566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm));
13399566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(linesearch->vec_func, NORM_2, &linesearch->fnorm));
13409566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
13419566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(linesearch->vec_update, NORM_2, &linesearch->ynorm));
13426a388c36SPeter Brune     }
13439bd66eb0SPeter Brune   }
13446a388c36SPeter Brune   PetscFunctionReturn(0);
13456a388c36SPeter Brune }
13466a388c36SPeter Brune 
13476f263ca3SPeter Brune /*@
13486f263ca3SPeter Brune    SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search.
13496f263ca3SPeter Brune 
13506f263ca3SPeter Brune    Input Parameters:
135178bcb3b5SPeter Brune +  linesearch  - linesearch context
135278bcb3b5SPeter Brune -  flg  - indicates whether or not to compute norms
13536f263ca3SPeter Brune 
13546f263ca3SPeter Brune    Options Database Keys:
13550b00b554SBarry Smith .   -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic (none) linesearch
13566f263ca3SPeter Brune 
1357*f6dfbefdSBarry Smith    Note:
1358*f6dfbefdSBarry 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.
13596f263ca3SPeter Brune 
13606f263ca3SPeter Brune    Level: intermediate
13616f263ca3SPeter Brune 
1362*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetNorms()`, `SNESLineSearchComputeNorms()`, `SNESLINESEARCHBASIC`
13636f263ca3SPeter Brune @*/
13649371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg) {
13656f263ca3SPeter Brune   PetscFunctionBegin;
13666f263ca3SPeter Brune   linesearch->norms = flg;
13676f263ca3SPeter Brune   PetscFunctionReturn(0);
13686f263ca3SPeter Brune }
13696f263ca3SPeter Brune 
1370f40b411bSPeter Brune /*@
1371*f6dfbefdSBarry Smith    SNESLineSearchGetVecs - Gets the vectors from the `SNESLineSearch` context
1372*f6dfbefdSBarry Smith 
1373*f6dfbefdSBarry Smith    Not Collective on but the vectors are parallel
1374f40b411bSPeter Brune 
1375f899ff85SJose E. Roman    Input Parameter:
137678bcb3b5SPeter Brune .  linesearch - linesearch context
1377f40b411bSPeter Brune 
1378f40b411bSPeter Brune    Output Parameters:
13796232e825SPeter Brune +  X - Solution vector
13806232e825SPeter Brune .  F - Function vector
13816232e825SPeter Brune .  Y - Search direction vector
13826232e825SPeter Brune .  W - Solution work vector
13836232e825SPeter Brune -  G - Function work vector
13846232e825SPeter Brune 
13857bba9028SPeter Brune    Notes:
13867bba9028SPeter Brune    At the beginning of a line search application, X should contain a
13876232e825SPeter Brune    solution and the vector F the function computed at X.  At the end of the
13886232e825SPeter Brune    line search application, X should contain the new solution, and F the
13896232e825SPeter Brune    function evaluated at the new solution.
1390f40b411bSPeter Brune 
1391*f6dfbefdSBarry Smith    These vectors are owned by the `SNESLineSearch` and should not be destroyed by the caller
13922a7a6963SBarry Smith 
139378bcb3b5SPeter Brune    Level: advanced
1394f40b411bSPeter Brune 
1395*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()`
1396f40b411bSPeter Brune @*/
13979371c9d4SSatish Balay PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch, Vec *X, Vec *F, Vec *Y, Vec *W, Vec *G) {
13986a388c36SPeter Brune   PetscFunctionBegin;
1399f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
14006a388c36SPeter Brune   if (X) {
14016a388c36SPeter Brune     PetscValidPointer(X, 2);
14026a388c36SPeter Brune     *X = linesearch->vec_sol;
14036a388c36SPeter Brune   }
14046a388c36SPeter Brune   if (F) {
14056a388c36SPeter Brune     PetscValidPointer(F, 3);
14066a388c36SPeter Brune     *F = linesearch->vec_func;
14076a388c36SPeter Brune   }
14086a388c36SPeter Brune   if (Y) {
14096a388c36SPeter Brune     PetscValidPointer(Y, 4);
14106a388c36SPeter Brune     *Y = linesearch->vec_update;
14116a388c36SPeter Brune   }
14126a388c36SPeter Brune   if (W) {
14136a388c36SPeter Brune     PetscValidPointer(W, 5);
14146a388c36SPeter Brune     *W = linesearch->vec_sol_new;
14156a388c36SPeter Brune   }
14166a388c36SPeter Brune   if (G) {
14176a388c36SPeter Brune     PetscValidPointer(G, 6);
14186a388c36SPeter Brune     *G = linesearch->vec_func_new;
14196a388c36SPeter Brune   }
14206a388c36SPeter Brune   PetscFunctionReturn(0);
14216a388c36SPeter Brune }
14226a388c36SPeter Brune 
1423f40b411bSPeter Brune /*@
1424*f6dfbefdSBarry Smith    SNESLineSearchSetVecs - Sets the vectors on the `SNESLineSearch` context
1425*f6dfbefdSBarry Smith 
1426*f6dfbefdSBarry Smith    Logically Collective on linesearch
1427f40b411bSPeter Brune 
1428f40b411bSPeter Brune    Input Parameters:
142978bcb3b5SPeter Brune +  linesearch - linesearch context
14306232e825SPeter Brune .  X - Solution vector
14316232e825SPeter Brune .  F - Function vector
14326232e825SPeter Brune .  Y - Search direction vector
14336232e825SPeter Brune .  W - Solution work vector
14346232e825SPeter Brune -  G - Function work vector
1435f40b411bSPeter Brune 
143678bcb3b5SPeter Brune    Level: advanced
1437f40b411bSPeter Brune 
1438*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetNorms()`, `SNESLineSearchGetVecs()`
1439f40b411bSPeter Brune @*/
14409371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch, Vec X, Vec F, Vec Y, Vec W, Vec G) {
14416a388c36SPeter Brune   PetscFunctionBegin;
1442f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
14436a388c36SPeter Brune   if (X) {
14446a388c36SPeter Brune     PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
14456a388c36SPeter Brune     linesearch->vec_sol = X;
14466a388c36SPeter Brune   }
14476a388c36SPeter Brune   if (F) {
14486a388c36SPeter Brune     PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
14496a388c36SPeter Brune     linesearch->vec_func = F;
14506a388c36SPeter Brune   }
14516a388c36SPeter Brune   if (Y) {
14526a388c36SPeter Brune     PetscValidHeaderSpecific(Y, VEC_CLASSID, 4);
14536a388c36SPeter Brune     linesearch->vec_update = Y;
14546a388c36SPeter Brune   }
14556a388c36SPeter Brune   if (W) {
14566a388c36SPeter Brune     PetscValidHeaderSpecific(W, VEC_CLASSID, 5);
14576a388c36SPeter Brune     linesearch->vec_sol_new = W;
14586a388c36SPeter Brune   }
14596a388c36SPeter Brune   if (G) {
14606a388c36SPeter Brune     PetscValidHeaderSpecific(G, VEC_CLASSID, 6);
14616a388c36SPeter Brune     linesearch->vec_func_new = G;
14626a388c36SPeter Brune   }
14636a388c36SPeter Brune   PetscFunctionReturn(0);
14646a388c36SPeter Brune }
14656a388c36SPeter Brune 
1466e7058c64SPeter Brune /*@C
1467f1c6b773SPeter Brune    SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1468*f6dfbefdSBarry Smith    `SNESLineSearch` options in the database.
1469e7058c64SPeter Brune 
1470*f6dfbefdSBarry Smith    Logically Collective on linesearch
1471e7058c64SPeter Brune 
1472e7058c64SPeter Brune    Input Parameters:
1473*f6dfbefdSBarry Smith +  linesearch - the `SNESLineSearch` context
1474e7058c64SPeter Brune -  prefix - the prefix to prepend to all option names
1475e7058c64SPeter Brune 
1476*f6dfbefdSBarry Smith    Note:
1477e7058c64SPeter Brune    A hyphen (-) must NOT be given at the beginning of the prefix name.
1478e7058c64SPeter Brune    The first character of all runtime options is AUTOMATICALLY the hyphen.
1479e7058c64SPeter Brune 
1480e7058c64SPeter Brune    Level: advanced
1481e7058c64SPeter Brune 
1482*f6dfbefdSBarry Smith .seealso: `SNESLineSearch()`, `SNESLineSearchSetFromOptions()`, `SNESGetOptionsPrefix()`
1483e7058c64SPeter Brune @*/
14849371c9d4SSatish Balay PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch, const char prefix[]) {
1485e7058c64SPeter Brune   PetscFunctionBegin;
1486f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
14879566063dSJacob Faibussowitsch   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)linesearch, prefix));
1488e7058c64SPeter Brune   PetscFunctionReturn(0);
1489e7058c64SPeter Brune }
1490e7058c64SPeter Brune 
1491e7058c64SPeter Brune /*@C
1492*f6dfbefdSBarry Smith    SNESLineSearchGetOptionsPrefix - Gets the prefix used for searching for all
1493f1c6b773SPeter Brune    SNESLineSearch options in the database.
1494e7058c64SPeter Brune 
1495e7058c64SPeter Brune    Not Collective
1496e7058c64SPeter Brune 
1497e7058c64SPeter Brune    Input Parameter:
1498*f6dfbefdSBarry Smith .  linesearch - the `SNESLineSearch` context
1499e7058c64SPeter Brune 
1500e7058c64SPeter Brune    Output Parameter:
1501e7058c64SPeter Brune .  prefix - pointer to the prefix string used
1502e7058c64SPeter Brune 
1503*f6dfbefdSBarry Smith    Fortran Note:
15048e557f58SPeter Brune    On the fortran side, the user should pass in a string 'prefix' of
1505e7058c64SPeter Brune    sufficient length to hold the prefix.
1506e7058c64SPeter Brune 
1507e7058c64SPeter Brune    Level: advanced
1508e7058c64SPeter Brune 
1509*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESAppendOptionsPrefix()`
1510e7058c64SPeter Brune @*/
15119371c9d4SSatish Balay PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch, const char *prefix[]) {
1512e7058c64SPeter Brune   PetscFunctionBegin;
1513f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
15149566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)linesearch, prefix));
1515e7058c64SPeter Brune   PetscFunctionReturn(0);
1516e7058c64SPeter Brune }
1517bf7f4e0aSPeter Brune 
15188d359177SBarry Smith /*@C
1519*f6dfbefdSBarry Smith    SNESLineSearchSetWorkVecs - Sets work vectors for the line search.
1520f40b411bSPeter Brune 
1521d8d19677SJose E. Roman    Input Parameters:
1522*f6dfbefdSBarry Smith +  linesearch - the `SNESLineSearch` context
1523f40b411bSPeter Brune -  nwork - the number of work vectors
1524f40b411bSPeter Brune 
1525f40b411bSPeter Brune    Level: developer
1526f40b411bSPeter Brune 
1527*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESSetWorkVecs()`
1528f40b411bSPeter Brune @*/
15299371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork) {
1530bf7f4e0aSPeter Brune   PetscFunctionBegin;
1531bf7f4e0aSPeter Brune   if (linesearch->vec_sol) {
15329566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work));
15338d359177SBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
1534bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1535bf7f4e0aSPeter Brune }
1536bf7f4e0aSPeter Brune 
1537f40b411bSPeter Brune /*@
1538422a814eSBarry Smith    SNESLineSearchGetReason - Gets the success/failure status of the last line search application
1539f40b411bSPeter Brune 
1540f40b411bSPeter Brune    Input Parameters:
154178bcb3b5SPeter Brune .  linesearch - linesearch context
1542f40b411bSPeter Brune 
1543f40b411bSPeter Brune    Output Parameters:
1544422a814eSBarry Smith .  result - The success or failure status
1545f40b411bSPeter Brune 
1546*f6dfbefdSBarry Smith    Note:
1547*f6dfbefdSBarry Smith    This is typically called after `SNESLineSearchApply()` in order to determine if the line-search failed
1548cd7522eaSPeter Brune    (and set the SNES convergence accordingly).
1549cd7522eaSPeter Brune 
1550*f6dfbefdSBarry Smith    Level: developer
1551f40b411bSPeter Brune 
1552*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetReason()`, `SNESLineSearchReason`
1553f40b411bSPeter Brune @*/
15549371c9d4SSatish Balay PetscErrorCode SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result) {
1555bf7f4e0aSPeter Brune   PetscFunctionBegin;
1556f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1557422a814eSBarry Smith   PetscValidPointer(result, 2);
1558422a814eSBarry Smith   *result = linesearch->result;
1559bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1560bf7f4e0aSPeter Brune }
1561bf7f4e0aSPeter Brune 
1562f40b411bSPeter Brune /*@
1563422a814eSBarry Smith    SNESLineSearchSetReason - Sets the success/failure status of the last line search application
1564f40b411bSPeter Brune 
1565f40b411bSPeter Brune    Input Parameters:
156678bcb3b5SPeter Brune +  linesearch - linesearch context
1567422a814eSBarry Smith -  result - The success or failure status
1568f40b411bSPeter Brune 
1569*f6dfbefdSBarry Smith    Note:
1570*f6dfbefdSBarry Smith    This is typically called in a `SNESLineSearchApply()` or a `SNESLINESEARCHSHELL` implementation to set
1571cd7522eaSPeter Brune    the success or failure of the line search method.
1572cd7522eaSPeter Brune 
157378bcb3b5SPeter Brune    Level: developer
1574f40b411bSPeter Brune 
1575*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchReason`, `SNESLineSearchGetSResult()`
1576f40b411bSPeter Brune @*/
15779371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result) {
15786a388c36SPeter Brune   PetscFunctionBegin;
15795fd66863SKarl Rupp   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1580422a814eSBarry Smith   linesearch->result = result;
15816a388c36SPeter Brune   PetscFunctionReturn(0);
15826a388c36SPeter Brune }
15836a388c36SPeter Brune 
15849bd66eb0SPeter Brune /*@C
1585f1c6b773SPeter Brune    SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.
15869bd66eb0SPeter Brune 
1587*f6dfbefdSBarry Smith    Logically Collective on linesearch
1588*f6dfbefdSBarry Smith 
15899bd66eb0SPeter Brune    Input Parameters:
1590*f6dfbefdSBarry Smith +  snes - nonlinear context obtained from `SNESCreate()`
15919bd66eb0SPeter Brune .  projectfunc - function for projecting the function to the bounds
15929bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
15939bd66eb0SPeter Brune 
15949bd66eb0SPeter Brune    Calling sequence of projectfunc:
15959bd66eb0SPeter Brune .vb
15969bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X)
15979bd66eb0SPeter Brune .ve
15989bd66eb0SPeter Brune 
15999bd66eb0SPeter Brune     Input parameters for projectfunc:
16009bd66eb0SPeter Brune +   snes - nonlinear context
16019bd66eb0SPeter Brune -   X - current solution
16029bd66eb0SPeter Brune 
1603*f6dfbefdSBarry Smith     Output parameter for projectfunc:
16049bd66eb0SPeter Brune .   X - Projected solution
16059bd66eb0SPeter Brune 
16069bd66eb0SPeter Brune    Calling sequence of normfunc:
16079bd66eb0SPeter Brune .vb
16089bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm)
16099bd66eb0SPeter Brune .ve
16109bd66eb0SPeter Brune 
1611cd7522eaSPeter Brune     Input parameters for normfunc:
16129bd66eb0SPeter Brune +   snes - nonlinear context
16139bd66eb0SPeter Brune .   X - current solution
16149bd66eb0SPeter Brune -   F - current residual
16159bd66eb0SPeter Brune 
1616*f6dfbefdSBarry Smith     Output parameter for normfunc:
16179bd66eb0SPeter Brune .   fnorm - VI-specific norm of the function
16189bd66eb0SPeter Brune 
1619*f6dfbefdSBarry Smith     Note:
1620cd7522eaSPeter Brune     The VI solvers require projection of the solution to the feasible set.  projectfunc should implement this.
1621cd7522eaSPeter Brune 
1622cd7522eaSPeter Brune     The VI solvers require special evaluation of the function norm such that the norm is only calculated
1623cd7522eaSPeter Brune     on the inactive set.  This should be implemented by normfunc.
16249bd66eb0SPeter Brune 
1625*f6dfbefdSBarry Smith     Level: advanced
16269bd66eb0SPeter Brune 
1627*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetVIFunctions()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()`
16289bd66eb0SPeter Brune @*/
16299371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc) {
16309bd66eb0SPeter Brune   PetscFunctionBegin;
1631f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
16329bd66eb0SPeter Brune   if (projectfunc) linesearch->ops->viproject = projectfunc;
16339bd66eb0SPeter Brune   if (normfunc) linesearch->ops->vinorm = normfunc;
16349bd66eb0SPeter Brune   PetscFunctionReturn(0);
16359bd66eb0SPeter Brune }
16369bd66eb0SPeter Brune 
16379bd66eb0SPeter Brune /*@C
1638f1c6b773SPeter Brune    SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.
16399bd66eb0SPeter Brune 
1640*f6dfbefdSBarry Smith    Not Collective
1641*f6dfbefdSBarry Smith 
1642f899ff85SJose E. Roman    Input Parameter:
1643*f6dfbefdSBarry Smith .  linesearch - the line search context, obtain with `SNESGetLineSearch()`
16449bd66eb0SPeter Brune 
16459bd66eb0SPeter Brune    Output Parameters:
16469bd66eb0SPeter Brune +  projectfunc - function for projecting the function to the bounds
16479bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
16489bd66eb0SPeter Brune 
1649*f6dfbefdSBarry Smith     Level: advanced
16509bd66eb0SPeter Brune 
1651*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetVIFunctions()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()`
16529bd66eb0SPeter Brune @*/
16539371c9d4SSatish Balay PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc) {
16549bd66eb0SPeter Brune   PetscFunctionBegin;
16559bd66eb0SPeter Brune   if (projectfunc) *projectfunc = linesearch->ops->viproject;
16569bd66eb0SPeter Brune   if (normfunc) *normfunc = linesearch->ops->vinorm;
16579bd66eb0SPeter Brune   PetscFunctionReturn(0);
16589bd66eb0SPeter Brune }
16599bd66eb0SPeter Brune 
1660bf7f4e0aSPeter Brune /*@C
1661*f6dfbefdSBarry Smith   SNESLineSearchRegister - register a line search method
1662bf7f4e0aSPeter Brune 
1663bf7f4e0aSPeter Brune   Level: advanced
1664*f6dfbefdSBarry Smith 
1665*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchSetType()`
1666bf7f4e0aSPeter Brune @*/
16679371c9d4SSatish Balay PetscErrorCode SNESLineSearchRegister(const char sname[], PetscErrorCode (*function)(SNESLineSearch)) {
1668bf7f4e0aSPeter Brune   PetscFunctionBegin;
16699566063dSJacob Faibussowitsch   PetscCall(SNESInitializePackage());
16709566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&SNESLineSearchList, sname, function));
1671bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1672bf7f4e0aSPeter Brune }
1673