xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision e4094ef18e7e53fda86cf35f3a47fda48a8e77d8)
1af0996ceSBarry Smith #include <petsc/private/linesearchimpl.h> /*I "petscsnes.h" I*/
2bf7f4e0aSPeter Brune 
3f1c6b773SPeter Brune PetscBool         SNESLineSearchRegisterAllCalled = PETSC_FALSE;
40298fd71SBarry Smith PetscFunctionList SNESLineSearchList              = NULL;
5bf7f4e0aSPeter Brune 
6f1c6b773SPeter Brune PetscClassId  SNESLINESEARCH_CLASSID;
757a83faaSBarry Smith PetscLogEvent SNESLINESEARCH_Apply;
8bf7f4e0aSPeter Brune 
9dcf2fd19SBarry Smith /*@
10f6dfbefdSBarry Smith   SNESLineSearchMonitorCancel - Clears all the monitor functions for a `SNESLineSearch` object.
11dcf2fd19SBarry Smith 
12c3339decSBarry Smith   Logically Collective
13dcf2fd19SBarry Smith 
142fe279fdSBarry Smith   Input Parameter:
15f6dfbefdSBarry Smith . ls - the `SNESLineSearch` context
16dcf2fd19SBarry Smith 
17dcf2fd19SBarry Smith   Options Database Key:
18dcf2fd19SBarry Smith . -snes_linesearch_monitor_cancel - cancels all monitors that have been hardwired
19f6dfbefdSBarry Smith     into a code by calls to `SNESLineSearchMonitorSet()`, but does not cancel those
20dcf2fd19SBarry Smith     set via the options database
21dcf2fd19SBarry Smith 
2220f4b53cSBarry Smith   Level: advanced
2320f4b53cSBarry Smith 
24dcf2fd19SBarry Smith   Notes:
25f6dfbefdSBarry Smith   There is no way to clear one specific monitor from a `SNESLineSearch` object.
26dcf2fd19SBarry Smith 
27f6dfbefdSBarry Smith   This does not clear the monitor set with `SNESLineSearchSetDefaultMonitor()` use `SNESLineSearchSetDefaultMonitor`(ls,NULL) to cancel
28dcf2fd19SBarry Smith   that one.
29dcf2fd19SBarry Smith 
30f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorDefault()`, `SNESLineSearchMonitorSet()`
31dcf2fd19SBarry Smith @*/
32d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchMonitorCancel(SNESLineSearch ls)
33d71ae5a4SJacob Faibussowitsch {
34dcf2fd19SBarry Smith   PetscInt i;
35dcf2fd19SBarry Smith 
36dcf2fd19SBarry Smith   PetscFunctionBegin;
37dcf2fd19SBarry Smith   PetscValidHeaderSpecific(ls, SNESLINESEARCH_CLASSID, 1);
38dcf2fd19SBarry Smith   for (i = 0; i < ls->numbermonitors; i++) {
3948a46eb9SPierre Jolivet     if (ls->monitordestroy[i]) PetscCall((*ls->monitordestroy[i])(&ls->monitorcontext[i]));
40dcf2fd19SBarry Smith   }
41dcf2fd19SBarry Smith   ls->numbermonitors = 0;
423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
43dcf2fd19SBarry Smith }
44dcf2fd19SBarry Smith 
45dcf2fd19SBarry Smith /*@
46dcf2fd19SBarry Smith   SNESLineSearchMonitor - runs the user provided monitor routines, if they exist
47dcf2fd19SBarry Smith 
48c3339decSBarry Smith   Collective
49dcf2fd19SBarry Smith 
502fe279fdSBarry Smith   Input Parameter:
51dcf2fd19SBarry Smith . ls - the linesearch object
52dcf2fd19SBarry Smith 
532fe279fdSBarry Smith   Level: developer
542fe279fdSBarry Smith 
55f6dfbefdSBarry Smith   Note:
5620f4b53cSBarry Smith   This routine is called by the `SNES` implementations.
57dcf2fd19SBarry Smith   It does not typically need to be called by the user.
58dcf2fd19SBarry Smith 
59f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorSet()`
60dcf2fd19SBarry Smith @*/
61d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchMonitor(SNESLineSearch ls)
62d71ae5a4SJacob Faibussowitsch {
63dcf2fd19SBarry Smith   PetscInt i, n = ls->numbermonitors;
64dcf2fd19SBarry Smith 
65dcf2fd19SBarry Smith   PetscFunctionBegin;
6648a46eb9SPierre Jolivet   for (i = 0; i < n; i++) PetscCall((*ls->monitorftns[i])(ls, ls->monitorcontext[i]));
673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
68dcf2fd19SBarry Smith }
69dcf2fd19SBarry Smith 
70dcf2fd19SBarry Smith /*@C
71dcf2fd19SBarry Smith   SNESLineSearchMonitorSet - Sets an ADDITIONAL function that is to be used at every
72dcf2fd19SBarry Smith   iteration of the nonlinear solver to display the iteration's
73dcf2fd19SBarry Smith   progress.
74dcf2fd19SBarry Smith 
75c3339decSBarry Smith   Logically Collective
76dcf2fd19SBarry Smith 
77dcf2fd19SBarry Smith   Input Parameters:
78f6dfbefdSBarry Smith + ls             - the `SNESLineSearch` context
79dcf2fd19SBarry Smith . f              - the monitor function
80dcf2fd19SBarry Smith . mctx           - [optional] user-defined context for private data for the
8120f4b53cSBarry Smith           monitor routine (use `NULL` if no context is desired)
82dcf2fd19SBarry Smith - monitordestroy - [optional] routine that frees monitor context
8320f4b53cSBarry Smith           (may be `NULL`)
8420f4b53cSBarry Smith 
8520f4b53cSBarry Smith   Level: intermediate
86dcf2fd19SBarry Smith 
87f6dfbefdSBarry Smith   Note:
88dcf2fd19SBarry Smith   Several different monitoring routines may be set by calling
89f6dfbefdSBarry Smith   `SNESLineSearchMonitorSet()` multiple times; all will be called in the
90dcf2fd19SBarry Smith   order in which they were set.
91dcf2fd19SBarry Smith 
92*e4094ef1SJacob Faibussowitsch   Fortran Notes:
93f6dfbefdSBarry Smith   Only a single monitor function can be set for each `SNESLineSearch` object
94dcf2fd19SBarry Smith 
95f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorDefault()`, `SNESLineSearchMonitorCancel()`
96dcf2fd19SBarry Smith @*/
97d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchMonitorSet(SNESLineSearch ls, PetscErrorCode (*f)(SNESLineSearch, void *), void *mctx, PetscErrorCode (*monitordestroy)(void **))
98d71ae5a4SJacob Faibussowitsch {
9978064530SBarry Smith   PetscInt  i;
10078064530SBarry Smith   PetscBool identical;
10178064530SBarry Smith 
102dcf2fd19SBarry Smith   PetscFunctionBegin;
103dcf2fd19SBarry Smith   PetscValidHeaderSpecific(ls, SNESLINESEARCH_CLASSID, 1);
10478064530SBarry Smith   for (i = 0; i < ls->numbermonitors; i++) {
1059566063dSJacob Faibussowitsch     PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))f, mctx, monitordestroy, (PetscErrorCode(*)(void))ls->monitorftns[i], ls->monitorcontext[i], ls->monitordestroy[i], &identical));
1063ba16761SJacob Faibussowitsch     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
10778064530SBarry Smith   }
1085f80ce2aSJacob Faibussowitsch   PetscCheck(ls->numbermonitors < MAXSNESLSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many monitors set");
109dcf2fd19SBarry Smith   ls->monitorftns[ls->numbermonitors]      = f;
110dcf2fd19SBarry Smith   ls->monitordestroy[ls->numbermonitors]   = monitordestroy;
111dcf2fd19SBarry Smith   ls->monitorcontext[ls->numbermonitors++] = (void *)mctx;
1123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
113dcf2fd19SBarry Smith }
114dcf2fd19SBarry Smith 
115dcf2fd19SBarry Smith /*@C
116f6dfbefdSBarry Smith   SNESLineSearchMonitorSolutionUpdate - Monitors each update of the function value the linesearch tries
117dcf2fd19SBarry Smith 
118c3339decSBarry Smith   Collective
119dcf2fd19SBarry Smith 
120dcf2fd19SBarry Smith   Input Parameters:
121f6dfbefdSBarry Smith + ls - the `SNES` linesearch object
122f6dfbefdSBarry Smith - vf - the context for the monitor, in this case it is an `PetscViewerAndFormat`
123dcf2fd19SBarry Smith 
124dcf2fd19SBarry Smith   Level: intermediate
125dcf2fd19SBarry Smith 
126f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESMonitorSet()`, `SNESMonitorSolution()`
127dcf2fd19SBarry Smith @*/
128d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchMonitorSolutionUpdate(SNESLineSearch ls, PetscViewerAndFormat *vf)
129d71ae5a4SJacob Faibussowitsch {
130d12e167eSBarry Smith   PetscViewer viewer = vf->viewer;
131dcf2fd19SBarry Smith   Vec         Y, W, G;
132dcf2fd19SBarry Smith 
133dcf2fd19SBarry Smith   PetscFunctionBegin;
1349566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetVecs(ls, NULL, NULL, &Y, &W, &G));
1359566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(viewer, vf->format));
1369566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted update to solution \n"));
1379566063dSJacob Faibussowitsch   PetscCall(VecView(Y, viewer));
1389566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted new solution \n"));
1399566063dSJacob Faibussowitsch   PetscCall(VecView(W, viewer));
1409566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted updated function value\n"));
1419566063dSJacob Faibussowitsch   PetscCall(VecView(G, viewer));
1429566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(viewer));
1433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
144dcf2fd19SBarry Smith }
145dcf2fd19SBarry Smith 
146f40b411bSPeter Brune /*@
147cd7522eaSPeter Brune   SNESLineSearchCreate - Creates the line search context.
148f40b411bSPeter Brune 
149f6dfbefdSBarry Smith   Logically Collective
150f40b411bSPeter Brune 
1512fe279fdSBarry Smith   Input Parameter:
152f6dfbefdSBarry Smith . comm - MPI communicator for the line search (typically from the associated `SNES` context).
153f40b411bSPeter Brune 
1542fe279fdSBarry Smith   Output Parameter:
1558e557f58SPeter Brune . outlinesearch - the new linesearch context
156f40b411bSPeter Brune 
157162e0bf5SPeter Brune   Level: developer
158162e0bf5SPeter Brune 
159f6dfbefdSBarry Smith   Note:
160f6dfbefdSBarry Smith   The preferred calling sequence for users is to use `SNESGetLineSearch()` to acquire the `SNESLineSearch` instance
161f6dfbefdSBarry Smith   already associated with the `SNES`.
162f40b411bSPeter Brune 
163f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `LineSearchDestroy()`, `SNESGetLineSearch()`
164f40b411bSPeter Brune @*/
165d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch)
166d71ae5a4SJacob Faibussowitsch {
167f1c6b773SPeter Brune   SNESLineSearch linesearch;
168bf388a1fSBarry Smith 
169bf7f4e0aSPeter Brune   PetscFunctionBegin;
170ea5d4fccSPeter Brune   PetscValidPointer(outlinesearch, 2);
1719566063dSJacob Faibussowitsch   PetscCall(SNESInitializePackage());
1720298fd71SBarry Smith   *outlinesearch = NULL;
173f5af7f23SKarl Rupp 
1749566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(linesearch, SNESLINESEARCH_CLASSID, "SNESLineSearch", "Linesearch", "SNESLineSearch", comm, SNESLineSearchDestroy, SNESLineSearchView));
175bf7f4e0aSPeter Brune 
1760298fd71SBarry Smith   linesearch->vec_sol_new  = NULL;
1770298fd71SBarry Smith   linesearch->vec_func_new = NULL;
1780298fd71SBarry Smith   linesearch->vec_sol      = NULL;
1790298fd71SBarry Smith   linesearch->vec_func     = NULL;
1800298fd71SBarry Smith   linesearch->vec_update   = NULL;
1819bd66eb0SPeter Brune 
182bf7f4e0aSPeter Brune   linesearch->lambda       = 1.0;
183bf7f4e0aSPeter Brune   linesearch->fnorm        = 1.0;
184bf7f4e0aSPeter Brune   linesearch->ynorm        = 1.0;
185bf7f4e0aSPeter Brune   linesearch->xnorm        = 1.0;
186422a814eSBarry Smith   linesearch->result       = SNES_LINESEARCH_SUCCEEDED;
187bf7f4e0aSPeter Brune   linesearch->norms        = PETSC_TRUE;
188bf7f4e0aSPeter Brune   linesearch->keeplambda   = PETSC_FALSE;
189bf7f4e0aSPeter Brune   linesearch->damping      = 1.0;
190bf7f4e0aSPeter Brune   linesearch->maxstep      = 1e8;
191bf7f4e0aSPeter Brune   linesearch->steptol      = 1e-12;
192516fe3c3SPeter Brune   linesearch->rtol         = 1e-8;
193516fe3c3SPeter Brune   linesearch->atol         = 1e-15;
194516fe3c3SPeter Brune   linesearch->ltol         = 1e-8;
1950298fd71SBarry Smith   linesearch->precheckctx  = NULL;
1960298fd71SBarry Smith   linesearch->postcheckctx = NULL;
19759405d5eSPeter Brune   linesearch->max_its      = 1;
198bf7f4e0aSPeter Brune   linesearch->setupcalled  = PETSC_FALSE;
1993add74b1SBarry Smith   linesearch->monitor      = NULL;
200bf7f4e0aSPeter Brune   *outlinesearch           = linesearch;
2013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
202bf7f4e0aSPeter Brune }
203bf7f4e0aSPeter Brune 
204f40b411bSPeter Brune /*@
20578bcb3b5SPeter Brune   SNESLineSearchSetUp - Prepares the line search for being applied by allocating
20678bcb3b5SPeter Brune   any required vectors.
207f40b411bSPeter Brune 
208c3339decSBarry Smith   Collective
209f40b411bSPeter Brune 
2102fe279fdSBarry Smith   Input Parameter:
211f6dfbefdSBarry Smith . linesearch - The `SNESLineSearch` instance.
212f40b411bSPeter Brune 
21320f4b53cSBarry Smith   Level: advanced
21420f4b53cSBarry Smith 
215f6dfbefdSBarry Smith   Note:
216f6dfbefdSBarry Smith   For most cases, this needn't be called by users or outside of `SNESLineSearchApply()`.
217cd7522eaSPeter Brune   The only current case where this is called outside of this is for the VI
21878bcb3b5SPeter Brune   solvers, which modify the solution and work vectors before the first call
219f6dfbefdSBarry Smith   of `SNESLineSearchApply()`, requiring the `SNESLineSearch` work vectors to be
220cd7522eaSPeter Brune   allocated upfront.
221cd7522eaSPeter Brune 
222f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchReset()`
223f40b411bSPeter Brune @*/
224d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch)
225d71ae5a4SJacob Faibussowitsch {
226bf7f4e0aSPeter Brune   PetscFunctionBegin;
22748a46eb9SPierre Jolivet   if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC));
228bf7f4e0aSPeter Brune   if (!linesearch->setupcalled) {
22948a46eb9SPierre Jolivet     if (!linesearch->vec_sol_new) PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new));
23048a46eb9SPierre Jolivet     if (!linesearch->vec_func_new) PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new));
231dbbe0bcdSBarry Smith     if (linesearch->ops->setup) PetscUseTypeMethod(linesearch, setup);
2329566063dSJacob Faibussowitsch     if (!linesearch->ops->snesfunc) PetscCall(SNESLineSearchSetFunction(linesearch, SNESComputeFunction));
233bf7f4e0aSPeter Brune     linesearch->lambda      = linesearch->damping;
234bf7f4e0aSPeter Brune     linesearch->setupcalled = PETSC_TRUE;
235bf7f4e0aSPeter Brune   }
2363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
237bf7f4e0aSPeter Brune }
238bf7f4e0aSPeter Brune 
239f40b411bSPeter Brune /*@
240f6dfbefdSBarry Smith   SNESLineSearchReset - Undoes the `SNESLineSearchSetUp()` and deletes any `Vec`s or `Mat`s allocated by the line search.
241f40b411bSPeter Brune 
242c3339decSBarry Smith   Collective
243f40b411bSPeter Brune 
2442fe279fdSBarry Smith   Input Parameter:
245f6dfbefdSBarry Smith . linesearch - The `SNESLineSearch` instance.
246f40b411bSPeter Brune 
24720f4b53cSBarry Smith   Level: developer
24820f4b53cSBarry Smith 
249f6dfbefdSBarry Smith   Note:
250f6dfbefdSBarry Smith   Usually only called by `SNESReset()`
251f190f2fcSBarry Smith 
252f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchSetUp()`
253f40b411bSPeter Brune @*/
254d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch)
255d71ae5a4SJacob Faibussowitsch {
256bf7f4e0aSPeter Brune   PetscFunctionBegin;
257dbbe0bcdSBarry Smith   if (linesearch->ops->reset) PetscUseTypeMethod(linesearch, reset);
258f5af7f23SKarl Rupp 
2599566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&linesearch->vec_sol_new));
2609566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&linesearch->vec_func_new));
261bf7f4e0aSPeter Brune 
2629566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(linesearch->nwork, &linesearch->work));
263f5af7f23SKarl Rupp 
264bf7f4e0aSPeter Brune   linesearch->nwork       = 0;
265bf7f4e0aSPeter Brune   linesearch->setupcalled = PETSC_FALSE;
2663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
267bf7f4e0aSPeter Brune }
268bf7f4e0aSPeter Brune 
269ed07d7d7SPeter Brune /*@C
270f6dfbefdSBarry Smith   SNESLineSearchSetFunction - Sets the function evaluation used by the `SNES` line search
271f6dfbefdSBarry Smith   `
272*e4094ef1SJacob Faibussowitsch 
273ed07d7d7SPeter Brune   Input Parameters:
274*e4094ef1SJacob Faibussowitsch + linesearch - the `SNESLineSearch` context
275*e4094ef1SJacob Faibussowitsch - func       - function evaluation routine, this is usually the function provided with `SNESSetFunction()`
276ed07d7d7SPeter Brune 
277ed07d7d7SPeter Brune   Level: developer
278ed07d7d7SPeter Brune 
279f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESSetFunction()`
280ed07d7d7SPeter Brune @*/
281d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES, Vec, Vec))
282d71ae5a4SJacob Faibussowitsch {
283ed07d7d7SPeter Brune   PetscFunctionBegin;
284ed07d7d7SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
285ed07d7d7SPeter Brune   linesearch->ops->snesfunc = func;
2863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
287ed07d7d7SPeter Brune }
288ed07d7d7SPeter Brune 
28986d74e61SPeter Brune /*@C
290f190f2fcSBarry Smith   SNESLineSearchSetPreCheck - Sets a user function that is called after the initial search direction has been computed but
291df3898eeSBarry Smith   before the line search routine has been applied. Allows the user to adjust the result of (usually a linear solve) that
292f190f2fcSBarry Smith   determined the search direction.
29386d74e61SPeter Brune 
294c3339decSBarry Smith   Logically Collective
29586d74e61SPeter Brune 
29686d74e61SPeter Brune   Input Parameters:
297f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context
29820f4b53cSBarry Smith . func       - [optional] function evaluation routine,  for the calling sequence see `SNESLineSearchPreCheck()`
29920f4b53cSBarry Smith - ctx        - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
30086d74e61SPeter Brune 
30186d74e61SPeter Brune   Level: intermediate
30286d74e61SPeter Brune 
303f6dfbefdSBarry Smith   Note:
304f0b84518SBarry Smith   Use `SNESLineSearchSetPostCheck()` to change the step after the line search.
305f0b84518SBarry Smith   search is complete.
306f0b84518SBarry Smith 
307f0b84518SBarry Smith   Use `SNESVISetVariableBounds()` and `SNESVISetComputeVariableBounds()` to cause `SNES` to automatically control the ranges of variables allowed.
308f0b84518SBarry Smith 
309f0b84518SBarry Smith .seealso: `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()`,
310869ba2dcSBarry Smith           `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`
311f0b84518SBarry Smith 
31286d74e61SPeter Brune @*/
313d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch, Vec, Vec, PetscBool *, void *), void *ctx)
314d71ae5a4SJacob Faibussowitsch {
3159bd66eb0SPeter Brune   PetscFunctionBegin;
316f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
317f190f2fcSBarry Smith   if (func) linesearch->ops->precheck = func;
31886d74e61SPeter Brune   if (ctx) linesearch->precheckctx = ctx;
3193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
32086d74e61SPeter Brune }
32186d74e61SPeter Brune 
32286d74e61SPeter Brune /*@C
32353deb899SBarry Smith   SNESLineSearchGetPreCheck - Gets the pre-check function for the line search routine.
32486d74e61SPeter Brune 
325f899ff85SJose E. Roman   Input Parameter:
326f6dfbefdSBarry Smith . linesearch - the `SNESLineSearch` context
32786d74e61SPeter Brune 
32886d74e61SPeter Brune   Output Parameters:
32920f4b53cSBarry Smith + func - [optional] function evaluation routine,  for calling sequence see `SNESLineSearchPreCheck`
33020f4b53cSBarry Smith - ctx  - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
33186d74e61SPeter Brune 
33286d74e61SPeter Brune   Level: intermediate
33386d74e61SPeter Brune 
334*e4094ef1SJacob Faibussowitsch .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchSetPostCheck()`
33586d74e61SPeter Brune @*/
336d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch, Vec, Vec, PetscBool *, void *), void **ctx)
337d71ae5a4SJacob Faibussowitsch {
3389bd66eb0SPeter Brune   PetscFunctionBegin;
339f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
340f190f2fcSBarry Smith   if (func) *func = linesearch->ops->precheck;
34186d74e61SPeter Brune   if (ctx) *ctx = linesearch->precheckctx;
3423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
34386d74e61SPeter Brune }
34486d74e61SPeter Brune 
34586d74e61SPeter Brune /*@C
346f190f2fcSBarry Smith   SNESLineSearchSetPostCheck - Sets a user function that is called after the line search has been applied to determine the step
347f190f2fcSBarry Smith   direction and length. Allows the user a chance to change or override the decision of the line search routine
34886d74e61SPeter Brune 
34920f4b53cSBarry Smith   Logically Collective
35086d74e61SPeter Brune 
35186d74e61SPeter Brune   Input Parameters:
352f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context
35320f4b53cSBarry Smith . func       - [optional] function evaluation routine,   for the calling sequence see `SNESLineSearchPostCheck`
35420f4b53cSBarry Smith - ctx        - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
35586d74e61SPeter Brune 
35686d74e61SPeter Brune   Level: intermediate
35786d74e61SPeter Brune 
358f0b84518SBarry Smith   Notes:
359f0b84518SBarry Smith   Use `SNESLineSearchSetPreCheck()` to change the step before the line search.
360f0b84518SBarry Smith   search is complete.
361f0b84518SBarry Smith 
362f0b84518SBarry Smith   Use `SNESVISetVariableBounds()` and `SNESVISetComputeVariableBounds()` to cause `SNES` to automatically control the ranges of variables allowed.
363f0b84518SBarry Smith 
364f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchGetPostCheck()`,
365*e4094ef1SJacob Faibussowitsch           `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`
36686d74e61SPeter Brune @*/
367d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx)
368d71ae5a4SJacob Faibussowitsch {
36986d74e61SPeter Brune   PetscFunctionBegin;
370f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
371f190f2fcSBarry Smith   if (func) linesearch->ops->postcheck = func;
37286d74e61SPeter Brune   if (ctx) linesearch->postcheckctx = ctx;
3733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
37486d74e61SPeter Brune }
37586d74e61SPeter Brune 
37686d74e61SPeter Brune /*@C
377f1c6b773SPeter Brune   SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine.
37886d74e61SPeter Brune 
379f899ff85SJose E. Roman   Input Parameter:
380f6dfbefdSBarry Smith . linesearch - the `SNESLineSearch` context
38186d74e61SPeter Brune 
38286d74e61SPeter Brune   Output Parameters:
383f6dfbefdSBarry Smith + func - [optional] function evaluation routine, see for the calling sequence `SNESLineSearchPostCheck`
38420f4b53cSBarry Smith - ctx  - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
38586d74e61SPeter Brune 
38686d74e61SPeter Brune   Level: intermediate
38786d74e61SPeter Brune 
388f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`
38986d74e61SPeter Brune @*/
390d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx)
391d71ae5a4SJacob Faibussowitsch {
3929bd66eb0SPeter Brune   PetscFunctionBegin;
393f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
394f190f2fcSBarry Smith   if (func) *func = linesearch->ops->postcheck;
39586d74e61SPeter Brune   if (ctx) *ctx = linesearch->postcheckctx;
3963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
39786d74e61SPeter Brune }
39886d74e61SPeter Brune 
399f40b411bSPeter Brune /*@
400f1c6b773SPeter Brune   SNESLineSearchPreCheck - Prepares the line search for being applied.
401f40b411bSPeter Brune 
402c3339decSBarry Smith   Logically Collective
403f40b411bSPeter Brune 
404f40b411bSPeter Brune   Input Parameters:
4057b1df9c1SPeter Brune + linesearch - The linesearch instance.
4067b1df9c1SPeter Brune . X          - The current solution
4077b1df9c1SPeter Brune - Y          - The step direction
408f40b411bSPeter Brune 
4092fe279fdSBarry Smith   Output Parameter:
4108e557f58SPeter Brune . changed - Indicator that the precheck routine has changed anything
411f40b411bSPeter Brune 
412f0b84518SBarry Smith   Level: advanced
413f40b411bSPeter Brune 
414f6dfbefdSBarry Smith   Note:
415f6dfbefdSBarry Smith   This calls any function provided with `SNESLineSearchSetPreCheck()`
416f6dfbefdSBarry Smith 
417f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`,
418f0b84518SBarry Smith           `SNESLineSearchGetPostCheck()``
419f40b411bSPeter Brune @*/
420d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed)
421d71ae5a4SJacob Faibussowitsch {
422bf7f4e0aSPeter Brune   PetscFunctionBegin;
423bf7f4e0aSPeter Brune   *changed = PETSC_FALSE;
4246b2b7091SBarry Smith   if (linesearch->ops->precheck) {
425dbbe0bcdSBarry Smith     PetscUseTypeMethod(linesearch, precheck, X, Y, changed, linesearch->precheckctx);
42638bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch, *changed, 4);
427bf7f4e0aSPeter Brune   }
4283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
429bf7f4e0aSPeter Brune }
430bf7f4e0aSPeter Brune 
431f40b411bSPeter Brune /*@
432ef46b1a6SStefano Zampini   SNESLineSearchPostCheck - Hook to modify step direction or updated solution after a successful linesearch
433f40b411bSPeter Brune 
434c3339decSBarry Smith   Logically Collective
435f40b411bSPeter Brune 
436f40b411bSPeter Brune   Input Parameters:
4377b1df9c1SPeter Brune + linesearch - The linesearch context
4387b1df9c1SPeter Brune . X          - The last solution
4397b1df9c1SPeter Brune . Y          - The step direction
4407b1df9c1SPeter Brune - W          - The updated solution, W = X + lambda*Y for some lambda
441f40b411bSPeter Brune 
442f40b411bSPeter Brune   Output Parameters:
44378bcb3b5SPeter Brune + changed_Y - Indicator if the direction Y has been changed.
44478bcb3b5SPeter Brune - changed_W - Indicator if the new candidate solution W has been changed.
445f40b411bSPeter Brune 
44620f4b53cSBarry Smith   Level: developer
44720f4b53cSBarry Smith 
448f6dfbefdSBarry Smith   Note:
449f6dfbefdSBarry Smith   This calls any function provided with `SNESLineSearchSetPreCheck()`
450f6dfbefdSBarry Smith 
451db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPrecheck()`, `SNESLineSearchGetPrecheck()`
452f40b411bSPeter Brune @*/
453d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W)
454d71ae5a4SJacob Faibussowitsch {
455bf7f4e0aSPeter Brune   PetscFunctionBegin;
456bf7f4e0aSPeter Brune   *changed_Y = PETSC_FALSE;
457bf7f4e0aSPeter Brune   *changed_W = PETSC_FALSE;
4586b2b7091SBarry Smith   if (linesearch->ops->postcheck) {
459dbbe0bcdSBarry Smith     PetscUseTypeMethod(linesearch, postcheck, X, Y, W, changed_Y, changed_W, linesearch->postcheckctx);
46038bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch, *changed_Y, 5);
46138bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch, *changed_W, 6);
46286d74e61SPeter Brune   }
4633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
46486d74e61SPeter Brune }
46586d74e61SPeter Brune 
46686d74e61SPeter Brune /*@C
46786d74e61SPeter Brune   SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration
46886d74e61SPeter Brune 
469c3339decSBarry Smith   Logically Collective
47086d74e61SPeter Brune 
4714165533cSJose E. Roman   Input Parameters:
47286d74e61SPeter Brune + linesearch - linesearch context
47386d74e61SPeter Brune . X          - base state for this step
474907376e6SBarry Smith - ctx        - context for this function
47586d74e61SPeter Brune 
47697bb3fdcSJose E. Roman   Input/Output Parameter:
47797bb3fdcSJose E. Roman . Y - correction, possibly modified
47897bb3fdcSJose E. Roman 
47997bb3fdcSJose E. Roman   Output Parameter:
48097bb3fdcSJose E. Roman . changed - flag indicating that Y was modified
48186d74e61SPeter Brune 
48286d74e61SPeter Brune   Options Database Key:
483cd7522eaSPeter Brune + -snes_linesearch_precheck_picard       - activate this routine
484cd7522eaSPeter Brune - -snes_linesearch_precheck_picard_angle - angle
48586d74e61SPeter Brune 
48686d74e61SPeter Brune   Level: advanced
48786d74e61SPeter Brune 
48886d74e61SPeter Brune   Notes:
489f6dfbefdSBarry Smith   This function should be passed to `SNESLineSearchSetPreCheck()`
49086d74e61SPeter Brune 
49186d74e61SPeter Brune   The justification for this method involves the linear convergence of a Picard iteration
49286d74e61SPeter Brune   so the Picard linearization should be provided in place of the "Jacobian". This correction
49386d74e61SPeter Brune   is generally not useful when using a Newton linearization.
49486d74e61SPeter Brune 
495*e4094ef1SJacob Faibussowitsch   References:
496f6dfbefdSBarry Smith   . - * - Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology.
49786d74e61SPeter Brune 
498f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESSetPicard()`, `SNESGetLineSearch()`, `SNESLineSearchSetPreCheck()`
49986d74e61SPeter Brune @*/
500d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed, void *ctx)
501d71ae5a4SJacob Faibussowitsch {
50286d74e61SPeter Brune   PetscReal   angle = *(PetscReal *)linesearch->precheckctx;
50386d74e61SPeter Brune   Vec         Ylast;
50486d74e61SPeter Brune   PetscScalar dot;
50586d74e61SPeter Brune   PetscInt    iter;
50686d74e61SPeter Brune   PetscReal   ynorm, ylastnorm, theta, angle_radians;
50786d74e61SPeter Brune   SNES        snes;
50886d74e61SPeter Brune 
50986d74e61SPeter Brune   PetscFunctionBegin;
5109566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
5119566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)snes, "SNESLineSearchPreCheckPicard_Ylast", (PetscObject *)&Ylast));
51286d74e61SPeter Brune   if (!Ylast) {
5139566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(Y, &Ylast));
5149566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)snes, "SNESLineSearchPreCheckPicard_Ylast", (PetscObject)Ylast));
5159566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)Ylast));
51686d74e61SPeter Brune   }
5179566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &iter));
51886d74e61SPeter Brune   if (iter < 2) {
5199566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y, Ylast));
52086d74e61SPeter Brune     *changed = PETSC_FALSE;
5213ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
52286d74e61SPeter Brune   }
52386d74e61SPeter Brune 
5249566063dSJacob Faibussowitsch   PetscCall(VecDot(Y, Ylast, &dot));
5259566063dSJacob Faibussowitsch   PetscCall(VecNorm(Y, NORM_2, &ynorm));
5269566063dSJacob Faibussowitsch   PetscCall(VecNorm(Ylast, NORM_2, &ylastnorm));
52786d74e61SPeter Brune   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
528255453a1SBarry Smith   theta         = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm), -1.0, 1.0));
52986d74e61SPeter Brune   angle_radians = angle * PETSC_PI / 180.;
53086d74e61SPeter Brune   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
53186d74e61SPeter Brune     /* Modify the step Y */
53286d74e61SPeter Brune     PetscReal alpha, ydiffnorm;
5339566063dSJacob Faibussowitsch     PetscCall(VecAXPY(Ylast, -1.0, Y));
5349566063dSJacob Faibussowitsch     PetscCall(VecNorm(Ylast, NORM_2, &ydiffnorm));
535f85e2ce2SBarry Smith     alpha = (ydiffnorm > .001 * ylastnorm) ? ylastnorm / ydiffnorm : 1000.0;
5369566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y, Ylast));
5379566063dSJacob Faibussowitsch     PetscCall(VecScale(Y, alpha));
5389566063dSJacob 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));
539a47ec85fSBarry Smith     *changed = PETSC_TRUE;
54086d74e61SPeter Brune   } else {
5419566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Angle %14.12e degrees exceeds threshold %14.12e, no correction applied\n", (double)(theta * 180. / PETSC_PI), (double)angle));
5429566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y, Ylast));
54386d74e61SPeter Brune     *changed = PETSC_FALSE;
544bf7f4e0aSPeter Brune   }
5453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
546bf7f4e0aSPeter Brune }
547bf7f4e0aSPeter Brune 
548f40b411bSPeter Brune /*@
549cd7522eaSPeter Brune   SNESLineSearchApply - Computes the line-search update.
550f40b411bSPeter Brune 
551c3339decSBarry Smith   Collective
552f40b411bSPeter Brune 
553f40b411bSPeter Brune   Input Parameters:
5548e557f58SPeter Brune + linesearch - The linesearch context
5558e557f58SPeter Brune - Y          - The search direction
556f40b411bSPeter Brune 
5576b867d5aSJose E. Roman   Input/Output Parameters:
5586b867d5aSJose E. Roman + X     - The current solution, on output the new solution
5596b867d5aSJose E. Roman . F     - The current function, on output the new function
5606b867d5aSJose E. Roman - fnorm - The current norm, on output the new function norm
561f40b411bSPeter Brune 
562cd7522eaSPeter Brune   Options Database Keys:
5630b00b554SBarry Smith + -snes_linesearch_type                - basic (or equivalently none), bt, l2, cp, nleqerr, shell
564dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches
5651fe24845SBarry Smith . -snes_linesearch_damping             - The linesearch damping parameter, default is 1.0 (no damping)
5661fe24845SBarry Smith . -snes_linesearch_norms               - Turn on/off the linesearch norms computation (SNESLineSearchSetComputeNorms())
5673c7d6663SPeter Brune . -snes_linesearch_keeplambda          - Keep the previous search length as the initial guess
5683c7d6663SPeter Brune - -snes_linesearch_max_it              - The number of iterations for iterative line searches
569cd7522eaSPeter Brune 
570*e4094ef1SJacob Faibussowitsch   Level: intermediate
57120f4b53cSBarry Smith 
572cd7522eaSPeter Brune   Notes:
573f6dfbefdSBarry Smith   This is typically called from within a `SNESSolve()` implementation in order to
574f6dfbefdSBarry Smith   help with convergence of the nonlinear method.  Various `SNES` types use line searches
575cd7522eaSPeter Brune   in different ways, but the overarching theme is that a line search is used to determine
576cd7522eaSPeter Brune   an optimal damping parameter of a step at each iteration of the method.  Each
577f6dfbefdSBarry Smith   application of the line search may invoke `SNESComputeFunction()` several times, and
578cd7522eaSPeter Brune   therefore may be fairly expensive.
579cd7522eaSPeter Brune 
580f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchPreCheck()`, `SNESLineSearchPostCheck()`, `SNESSolve()`, `SNESComputeFunction()`, `SNESLineSearchSetComputeNorms()`,
581db781477SPatrick Sanan           `SNESLineSearchType`, `SNESLineSearchSetType()`
582f40b411bSPeter Brune @*/
583d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal *fnorm, Vec Y)
584d71ae5a4SJacob Faibussowitsch {
585bf388a1fSBarry Smith   PetscFunctionBegin;
586f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
587bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
588bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
589064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(Y, VEC_CLASSID, 5);
590bf7f4e0aSPeter Brune 
591422a814eSBarry Smith   linesearch->result = SNES_LINESEARCH_SUCCEEDED;
592bf7f4e0aSPeter Brune 
593bf7f4e0aSPeter Brune   linesearch->vec_sol    = X;
594bf7f4e0aSPeter Brune   linesearch->vec_update = Y;
595bf7f4e0aSPeter Brune   linesearch->vec_func   = F;
596bf7f4e0aSPeter Brune 
5979566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchSetUp(linesearch));
598bf7f4e0aSPeter Brune 
599f5af7f23SKarl Rupp   if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */
600bf7f4e0aSPeter Brune 
601f5af7f23SKarl Rupp   if (fnorm) linesearch->fnorm = *fnorm;
60248a46eb9SPierre Jolivet   else PetscCall(VecNorm(F, NORM_2, &linesearch->fnorm));
603bf7f4e0aSPeter Brune 
6049566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(SNESLINESEARCH_Apply, linesearch, X, F, Y));
605bf7f4e0aSPeter Brune 
606dbbe0bcdSBarry Smith   PetscUseTypeMethod(linesearch, apply);
607bf7f4e0aSPeter Brune 
6089566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(SNESLINESEARCH_Apply, linesearch, X, F, Y));
609bf7f4e0aSPeter Brune 
610f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
6113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
612bf7f4e0aSPeter Brune }
613bf7f4e0aSPeter Brune 
614f40b411bSPeter Brune /*@
615f1c6b773SPeter Brune   SNESLineSearchDestroy - Destroys the line search instance.
616f40b411bSPeter Brune 
617c3339decSBarry Smith   Collective
618f40b411bSPeter Brune 
6192fe279fdSBarry Smith   Input Parameter:
6208e557f58SPeter Brune . linesearch - The linesearch context
621f40b411bSPeter Brune 
62284238204SBarry Smith   Level: developer
623f40b411bSPeter Brune 
624f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchReset()`, `SNESDestroy()`
625f40b411bSPeter Brune @*/
626d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchDestroy(SNESLineSearch *linesearch)
627d71ae5a4SJacob Faibussowitsch {
628bf7f4e0aSPeter Brune   PetscFunctionBegin;
6293ba16761SJacob Faibussowitsch   if (!*linesearch) PetscFunctionReturn(PETSC_SUCCESS);
630f1c6b773SPeter Brune   PetscValidHeaderSpecific((*linesearch), SNESLINESEARCH_CLASSID, 1);
6319371c9d4SSatish Balay   if (--((PetscObject)(*linesearch))->refct > 0) {
6329371c9d4SSatish Balay     *linesearch = NULL;
6333ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
6349371c9d4SSatish Balay   }
6359566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsViewOff((PetscObject)*linesearch));
6369566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchReset(*linesearch));
6373ba16761SJacob Faibussowitsch   PetscTryTypeMethod(*linesearch, destroy);
6389566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&(*linesearch)->monitor));
6399566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchMonitorCancel((*linesearch)));
6409566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(linesearch));
6413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
642bf7f4e0aSPeter Brune }
643bf7f4e0aSPeter Brune 
644f40b411bSPeter Brune /*@
645dcf2fd19SBarry Smith   SNESLineSearchSetDefaultMonitor - Turns on/off printing useful information and debugging output about the line search.
646bf7f4e0aSPeter Brune 
647c3339decSBarry Smith   Logically Collective
648f6dfbefdSBarry Smith 
649bf7f4e0aSPeter Brune   Input Parameters:
650dcf2fd19SBarry Smith + linesearch - the linesearch object
65120f4b53cSBarry Smith - viewer     - an `PETSCVIEWERASCII` `PetscViewer` or `NULL` to turn off monitor
652bf7f4e0aSPeter Brune 
653f6dfbefdSBarry Smith   Options Database Key:
654dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - enables the monitor
655bf7f4e0aSPeter Brune 
656bf7f4e0aSPeter Brune   Level: intermediate
657bf7f4e0aSPeter Brune 
658*e4094ef1SJacob Faibussowitsch   Developer Notes:
659f6dfbefdSBarry Smith   This monitor is implemented differently than the other line search monitors that are set with
660f6dfbefdSBarry Smith   `SNESLineSearchMonitorSet()` since it is called in many locations of the line search routines to display aspects of the
661d12e167eSBarry Smith   line search that are not visible to the other monitors.
662dcf2fd19SBarry Smith 
663f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `PETSCVIEWERASCII`, `SNESGetLineSearch()`, `SNESLineSearchGetDefaultMonitor()`, `PetscViewer`, `SNESLineSearchSetMonitor()`,
664f6dfbefdSBarry Smith           `SNESLineSearchMonitorSetFromOptions()`
665bf7f4e0aSPeter Brune @*/
666d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetDefaultMonitor(SNESLineSearch linesearch, PetscViewer viewer)
667d71ae5a4SJacob Faibussowitsch {
668bf7f4e0aSPeter Brune   PetscFunctionBegin;
6699566063dSJacob Faibussowitsch   if (viewer) PetscCall(PetscObjectReference((PetscObject)viewer));
6709566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&linesearch->monitor));
671dcf2fd19SBarry Smith   linesearch->monitor = viewer;
6723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
673bf7f4e0aSPeter Brune }
674bf7f4e0aSPeter Brune 
675f40b411bSPeter Brune /*@
676f6dfbefdSBarry Smith   SNESLineSearchGetDefaultMonitor - Gets the `PetscViewer` instance for the line search monitor.
677f6dfbefdSBarry Smith 
678c3339decSBarry Smith   Logically Collective
6796a388c36SPeter Brune 
680f190f2fcSBarry Smith   Input Parameter:
6818e557f58SPeter Brune . linesearch - linesearch context
682f40b411bSPeter Brune 
683f190f2fcSBarry Smith   Output Parameter:
6848e557f58SPeter Brune . monitor - monitor context
685f40b411bSPeter Brune 
686f40b411bSPeter Brune   Level: intermediate
687f40b411bSPeter Brune 
688f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchSetDefaultMonitor()`, `PetscViewer`
689f40b411bSPeter Brune @*/
690d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetDefaultMonitor(SNESLineSearch linesearch, PetscViewer *monitor)
691d71ae5a4SJacob Faibussowitsch {
6926a388c36SPeter Brune   PetscFunctionBegin;
693f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
6946a388c36SPeter Brune   *monitor = linesearch->monitor;
6953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6966a388c36SPeter Brune }
6976a388c36SPeter Brune 
698dcf2fd19SBarry Smith /*@C
699dcf2fd19SBarry Smith   SNESLineSearchMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
700dcf2fd19SBarry Smith 
701c3339decSBarry Smith   Collective
702dcf2fd19SBarry Smith 
703dcf2fd19SBarry Smith   Input Parameters:
704f6dfbefdSBarry Smith + ls           - `SNESLineSearch` object you wish to monitor
705dcf2fd19SBarry Smith . name         - the monitor type one is seeking
706dcf2fd19SBarry Smith . help         - message indicating what monitoring is done
707dcf2fd19SBarry Smith . manual       - manual page for the monitor
708dcf2fd19SBarry Smith . monitor      - the monitor function
709f6dfbefdSBarry 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`
710dcf2fd19SBarry Smith 
711dcf2fd19SBarry Smith   Level: developer
712dcf2fd19SBarry Smith 
713f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetMonitor()`, `PetscOptionsGetViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
714db781477SPatrick Sanan           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
715*e4094ef1SJacob Faibussowitsch           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
716db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
717c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
718db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
719db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
720dcf2fd19SBarry Smith @*/
721d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchMonitorSetFromOptions(SNESLineSearch ls, const char name[], const char help[], const char manual[], PetscErrorCode (*monitor)(SNESLineSearch, PetscViewerAndFormat *), PetscErrorCode (*monitorsetup)(SNESLineSearch, PetscViewerAndFormat *))
722d71ae5a4SJacob Faibussowitsch {
723dcf2fd19SBarry Smith   PetscViewer       viewer;
724dcf2fd19SBarry Smith   PetscViewerFormat format;
725dcf2fd19SBarry Smith   PetscBool         flg;
726dcf2fd19SBarry Smith 
727dcf2fd19SBarry Smith   PetscFunctionBegin;
7289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)ls), ((PetscObject)ls)->options, ((PetscObject)ls)->prefix, name, &viewer, &format, &flg));
729dcf2fd19SBarry Smith   if (flg) {
730d12e167eSBarry Smith     PetscViewerAndFormat *vf;
7319566063dSJacob Faibussowitsch     PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf));
7329566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)viewer));
7331baa6e33SBarry Smith     if (monitorsetup) PetscCall((*monitorsetup)(ls, vf));
7349566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchMonitorSet(ls, (PetscErrorCode(*)(SNESLineSearch, void *))monitor, vf, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy));
735dcf2fd19SBarry Smith   }
7363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
737dcf2fd19SBarry Smith }
738dcf2fd19SBarry Smith 
739f40b411bSPeter Brune /*@
740f1c6b773SPeter Brune   SNESLineSearchSetFromOptions - Sets options for the line search
741f40b411bSPeter Brune 
742c3339decSBarry Smith   Logically Collective
743f6dfbefdSBarry Smith 
744f899ff85SJose E. Roman   Input Parameter:
7458e557f58SPeter Brune . linesearch - linesearch context
746f40b411bSPeter Brune 
747f40b411bSPeter Brune   Options Database Keys:
7480b00b554SBarry Smith + -snes_linesearch_type <type>                                      - basic (or equivalently none), bt, l2, cp, nleqerr, shell
7493c7d6663SPeter Brune . -snes_linesearch_order <order>                                    - 1, 2, 3.  Most types only support certain orders (bt supports 2 or 3)
750f6dfbefdSBarry Smith . -snes_linesearch_norms                                            - Turn on/off the linesearch norms for the basic linesearch typem (`SNESLineSearchSetComputeNorms()`)
75171eef1aeSPeter Brune . -snes_linesearch_minlambda                                        - The minimum step length
7521a9b3a06SPeter Brune . -snes_linesearch_maxstep                                          - The maximum step size
7531a9b3a06SPeter Brune . -snes_linesearch_rtol                                             - Relative tolerance for iterative line searches
7541a9b3a06SPeter Brune . -snes_linesearch_atol                                             - Absolute tolerance for iterative line searches
7551a9b3a06SPeter Brune . -snes_linesearch_ltol                                             - Change in lambda tolerance for iterative line searches
7561a9b3a06SPeter Brune . -snes_linesearch_max_it                                           - The number of iterations for iterative line searches
757dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename]                              - Print progress of line searches
758dcf2fd19SBarry Smith . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine
7598e557f58SPeter Brune . -snes_linesearch_damping                                          - The linesearch damping parameter
760cd7522eaSPeter Brune . -snes_linesearch_keeplambda                                       - Keep the previous search length as the initial guess.
7611a9b3a06SPeter Brune . -snes_linesearch_precheck_picard                                  - Use precheck that speeds up convergence of picard method
762d8d34be6SBarry Smith - -snes_linesearch_precheck_picard_angle                            - Angle used in Picard precheck method
763f40b411bSPeter Brune 
764f40b411bSPeter Brune   Level: intermediate
765f40b411bSPeter Brune 
766f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchCreate()`, `SNESLineSearchSetOrder()`, `SNESLineSearchSetType()`, `SNESLineSearchSetTolerances()`, `SNESLineSearchSetDamping()`, `SNESLineSearchPreCheckPicard()`,
767db781477SPatrick Sanan           `SNESLineSearchType`, `SNESLineSearchSetComputeNorms()`
768f40b411bSPeter Brune @*/
769d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch)
770d71ae5a4SJacob Faibussowitsch {
7711a4f838cSPeter Brune   const char *deft = SNESLINESEARCHBASIC;
772bf7f4e0aSPeter Brune   char        type[256];
773bf7f4e0aSPeter Brune   PetscBool   flg, set;
774dcf2fd19SBarry Smith   PetscViewer viewer;
775bf388a1fSBarry Smith 
776bf7f4e0aSPeter Brune   PetscFunctionBegin;
7779566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchRegisterAll());
778bf7f4e0aSPeter Brune 
779d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)linesearch);
780f5af7f23SKarl Rupp   if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name;
7819566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-snes_linesearch_type", "Linesearch type", "SNESLineSearchSetType", SNESLineSearchList, deft, type, 256, &flg));
782bf7f4e0aSPeter Brune   if (flg) {
7839566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetType(linesearch, type));
784bf7f4e0aSPeter Brune   } else if (!((PetscObject)linesearch)->type_name) {
7859566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetType(linesearch, deft));
786bf7f4e0aSPeter Brune   }
787bf7f4e0aSPeter Brune 
7889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)linesearch), ((PetscObject)linesearch)->options, ((PetscObject)linesearch)->prefix, "-snes_linesearch_monitor", &viewer, NULL, &set));
789dcf2fd19SBarry Smith   if (set) {
7909566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetDefaultMonitor(linesearch, viewer));
7919566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
792dcf2fd19SBarry Smith   }
7939566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchMonitorSetFromOptions(linesearch, "-snes_linesearch_monitor_solution_update", "View correction at each iteration", "SNESLineSearchMonitorSolutionUpdate", SNESLineSearchMonitorSolutionUpdate, NULL));
794bf7f4e0aSPeter Brune 
7951a9b3a06SPeter Brune   /* tolerances */
7969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_minlambda", "Minimum step length", "SNESLineSearchSetTolerances", linesearch->steptol, &linesearch->steptol, NULL));
7979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_maxstep", "Maximum step size", "SNESLineSearchSetTolerances", linesearch->maxstep, &linesearch->maxstep, NULL));
7989566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_rtol", "Relative tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->rtol, &linesearch->rtol, NULL));
7999566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_atol", "Absolute tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->atol, &linesearch->atol, NULL));
8009566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_ltol", "Change in lambda tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->ltol, &linesearch->ltol, NULL));
8019566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_linesearch_max_it", "Maximum iterations for iterative line searches", "SNESLineSearchSetTolerances", linesearch->max_its, &linesearch->max_its, NULL));
8027a35526eSPeter Brune 
8031a9b3a06SPeter Brune   /* damping parameters */
8049566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_damping", "Line search damping and initial step guess", "SNESLineSearchSetDamping", linesearch->damping, &linesearch->damping, NULL));
8051a9b3a06SPeter Brune 
8069566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_linesearch_keeplambda", "Use previous lambda as damping", "SNESLineSearchSetKeepLambda", linesearch->keeplambda, &linesearch->keeplambda, NULL));
8071a9b3a06SPeter Brune 
8081a9b3a06SPeter Brune   /* precheck */
8099566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_linesearch_precheck_picard", "Use a correction that sometimes improves convergence of Picard iteration", "SNESLineSearchPreCheckPicard", flg, &flg, &set));
8107a35526eSPeter Brune   if (set) {
8117a35526eSPeter Brune     if (flg) {
8127a35526eSPeter Brune       linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */
813f5af7f23SKarl Rupp 
814d0609cedSBarry 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));
8159566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchSetPreCheck(linesearch, SNESLineSearchPreCheckPicard, &linesearch->precheck_picard_angle));
8167a35526eSPeter Brune     } else {
8179566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchSetPreCheck(linesearch, NULL, NULL));
8187a35526eSPeter Brune     }
8197a35526eSPeter Brune   }
8209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_linesearch_order", "Order of approximation used in the line search", "SNESLineSearchSetOrder", linesearch->order, &linesearch->order, NULL));
8219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_linesearch_norms", "Compute final norms in line search", "SNESLineSearchSetComputeNorms", linesearch->norms, &linesearch->norms, NULL));
8227a35526eSPeter Brune 
823dbbe0bcdSBarry Smith   PetscTryTypeMethod(linesearch, setfromoptions, PetscOptionsObject);
8245a9b6599SPeter Brune 
825dbbe0bcdSBarry Smith   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)linesearch, PetscOptionsObject));
826d0609cedSBarry Smith   PetscOptionsEnd();
8273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
828bf7f4e0aSPeter Brune }
829bf7f4e0aSPeter Brune 
830f40b411bSPeter Brune /*@
831f190f2fcSBarry Smith   SNESLineSearchView - Prints useful information about the line search
832f40b411bSPeter Brune 
83320f4b53cSBarry Smith   Logically Collective
83420f4b53cSBarry Smith 
835f40b411bSPeter Brune   Input Parameters:
8362fe279fdSBarry Smith + linesearch - linesearch context
8372fe279fdSBarry Smith - viewer     - the viewer to display the line search information
838f40b411bSPeter Brune 
839f40b411bSPeter Brune   Level: intermediate
840f40b411bSPeter Brune 
841f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `PetscViewer`, `SNESLineSearchCreate()`
842f40b411bSPeter Brune @*/
843d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer)
844d71ae5a4SJacob Faibussowitsch {
8457f1410a3SPeter Brune   PetscBool iascii;
846bf388a1fSBarry Smith 
847bf7f4e0aSPeter Brune   PetscFunctionBegin;
8487f1410a3SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
84948a46eb9SPierre Jolivet   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch), &viewer));
8507f1410a3SPeter Brune   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
8517f1410a3SPeter Brune   PetscCheckSameComm(linesearch, 1, viewer, 2);
852f40b411bSPeter Brune 
8539566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
8547f1410a3SPeter Brune   if (iascii) {
8559566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)linesearch, viewer));
8569566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
857dbbe0bcdSBarry Smith     PetscTryTypeMethod(linesearch, view, viewer);
8589566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
8599566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep, (double)linesearch->steptol));
8609566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol, (double)linesearch->atol, (double)linesearch->ltol));
86163a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  maximum iterations=%" PetscInt_FMT "\n", linesearch->max_its));
8626b2b7091SBarry Smith     if (linesearch->ops->precheck) {
8636b2b7091SBarry Smith       if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) {
86463a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  using precheck step to speed up Picard convergence\n"));
8657f1410a3SPeter Brune       } else {
86663a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  using user-defined precheck step\n"));
8677f1410a3SPeter Brune       }
8687f1410a3SPeter Brune     }
86948a46eb9SPierre Jolivet     if (linesearch->ops->postcheck) PetscCall(PetscViewerASCIIPrintf(viewer, "  using user-defined postcheck step\n"));
8707f1410a3SPeter Brune   }
8713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
872bf7f4e0aSPeter Brune }
873bf7f4e0aSPeter Brune 
874ea5d4fccSPeter Brune /*@C
875a80ff896SJed Brown   SNESLineSearchGetType - Gets the linesearch type
876a80ff896SJed Brown 
877c3339decSBarry Smith   Logically Collective
878a80ff896SJed Brown 
8792fe279fdSBarry Smith   Input Parameter:
880a80ff896SJed Brown . linesearch - linesearch context
881a80ff896SJed Brown 
8822fe279fdSBarry Smith   Output Parameter:
8832fe279fdSBarry Smith . type - The type of line search, or `NULL` if not set
884a80ff896SJed Brown 
885a80ff896SJed Brown   Level: intermediate
886a80ff896SJed Brown 
887*e4094ef1SJacob Faibussowitsch .seealso: `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchSetType()`
888a80ff896SJed Brown @*/
889d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetType(SNESLineSearch linesearch, SNESLineSearchType *type)
890d71ae5a4SJacob Faibussowitsch {
891a80ff896SJed Brown   PetscFunctionBegin;
892a80ff896SJed Brown   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
893dadcf809SJacob Faibussowitsch   PetscValidPointer(type, 2);
894a80ff896SJed Brown   *type = ((PetscObject)linesearch)->type_name;
8953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
896a80ff896SJed Brown }
897a80ff896SJed Brown 
898a80ff896SJed Brown /*@C
899f1c6b773SPeter Brune   SNESLineSearchSetType - Sets the linesearch type
900f40b411bSPeter Brune 
901c3339decSBarry Smith   Logically Collective
902f190f2fcSBarry Smith 
903f40b411bSPeter Brune   Input Parameters:
9048e557f58SPeter Brune + linesearch - linesearch context
905f40b411bSPeter Brune - type       - The type of line search to be used
906f40b411bSPeter Brune 
907cd7522eaSPeter Brune   Available Types:
908f6dfbefdSBarry Smith +  `SNESLINESEARCHBASIC` - (or equivalently `SNESLINESEARCHNONE`) Simple damping line search, defaults to using the full Newton step
909f6dfbefdSBarry Smith .  `SNESLINESEARCHBT` - Backtracking line search over the L2 norm of the function
910f6dfbefdSBarry Smith .  `SNESLINESEARCHL2` - Secant line search over the L2 norm of the function
911f6dfbefdSBarry Smith .  `SNESLINESEARCHCP` - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x)
912f6dfbefdSBarry Smith .  `SNESLINESEARCHNLEQERR` - Affine-covariant error-oriented linesearch
913f6dfbefdSBarry Smith -  `SNESLINESEARCHSHELL` - User provided `SNESLineSearch` implementation
9141fe24845SBarry Smith 
9153c7db156SBarry Smith   Options Database Key:
9160b00b554SBarry Smith . -snes_linesearch_type <type> - basic (or equivalently none), bt, l2, cp, nleqerr, shell
917cd7522eaSPeter Brune 
918f40b411bSPeter Brune   Level: intermediate
919f40b411bSPeter Brune 
920*e4094ef1SJacob Faibussowitsch .seealso: `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchGetType()`
921f40b411bSPeter Brune @*/
922d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type)
923d71ae5a4SJacob Faibussowitsch {
924bf7f4e0aSPeter Brune   PetscBool match;
9255f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(SNESLineSearch);
926bf7f4e0aSPeter Brune 
927bf7f4e0aSPeter Brune   PetscFunctionBegin;
928f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
929bf7f4e0aSPeter Brune   PetscValidCharPointer(type, 2);
930bf7f4e0aSPeter Brune 
9319566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)linesearch, type, &match));
9323ba16761SJacob Faibussowitsch   if (match) PetscFunctionReturn(PETSC_SUCCESS);
933bf7f4e0aSPeter Brune 
9349566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(SNESLineSearchList, type, &r));
9356adde796SStefano Zampini   PetscCheck(r, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested Line Search type %s", type);
936bf7f4e0aSPeter Brune   /* Destroy the previous private linesearch context */
937bf7f4e0aSPeter Brune   if (linesearch->ops->destroy) {
9389566063dSJacob Faibussowitsch     PetscCall((*(linesearch)->ops->destroy)(linesearch));
9390298fd71SBarry Smith     linesearch->ops->destroy = NULL;
940bf7f4e0aSPeter Brune   }
941f1c6b773SPeter Brune   /* Reinitialize function pointers in SNESLineSearchOps structure */
9429e5d0892SLisandro Dalcin   linesearch->ops->apply          = NULL;
9439e5d0892SLisandro Dalcin   linesearch->ops->view           = NULL;
9449e5d0892SLisandro Dalcin   linesearch->ops->setfromoptions = NULL;
9459e5d0892SLisandro Dalcin   linesearch->ops->destroy        = NULL;
946bf7f4e0aSPeter Brune 
9479566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)linesearch, type));
9489566063dSJacob Faibussowitsch   PetscCall((*r)(linesearch));
9493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
950bf7f4e0aSPeter Brune }
951bf7f4e0aSPeter Brune 
952f40b411bSPeter Brune /*@
953f6dfbefdSBarry Smith   SNESLineSearchSetSNES - Sets the `SNES` for the linesearch for function evaluation.
954f40b411bSPeter Brune 
955f40b411bSPeter Brune   Input Parameters:
9568e557f58SPeter Brune + linesearch - linesearch context
957f40b411bSPeter Brune - snes       - The snes instance
958f40b411bSPeter Brune 
95978bcb3b5SPeter Brune   Level: developer
96078bcb3b5SPeter Brune 
961f6dfbefdSBarry Smith   Note:
962f190f2fcSBarry Smith   This happens automatically when the line search is obtained/created with
963f6dfbefdSBarry Smith   `SNESGetLineSearch()`.  This routine is therefore mainly called within `SNES`
96478bcb3b5SPeter Brune   implementations.
965f40b411bSPeter Brune 
966f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetSNES()`, `SNESLineSearchSetVecs()`, `SNES`
967f40b411bSPeter Brune @*/
968d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes)
969d71ae5a4SJacob Faibussowitsch {
970bf7f4e0aSPeter Brune   PetscFunctionBegin;
971f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
972bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 2);
973bf7f4e0aSPeter Brune   linesearch->snes = snes;
9743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
975bf7f4e0aSPeter Brune }
976bf7f4e0aSPeter Brune 
977f40b411bSPeter Brune /*@
978f6dfbefdSBarry Smith   SNESLineSearchGetSNES - Gets the `SNES` instance associated with the line search.
979f6dfbefdSBarry Smith   Having an associated `SNES` is necessary because most line search implementations must be able to
980f6dfbefdSBarry Smith   evaluate the function using `SNESComputeFunction()` for the associated `SNES`.  This routine
981f6dfbefdSBarry Smith   is used in the line search implementations when one must get this associated `SNES` instance.
982f6dfbefdSBarry Smith 
983f6dfbefdSBarry Smith   Not Collective
984f40b411bSPeter Brune 
9852fe279fdSBarry Smith   Input Parameter:
9868e557f58SPeter Brune . linesearch - linesearch context
987f40b411bSPeter Brune 
9882fe279fdSBarry Smith   Output Parameter:
9892fe279fdSBarry Smith . snes - The `SNES` instance
990f40b411bSPeter Brune 
9918141a3b9SPeter Brune   Level: developer
992f40b411bSPeter Brune 
993f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESType`, `SNESLineSearchGetSNES()`, `SNESLineSearchSetVecs()`, `SNES`
994f40b411bSPeter Brune @*/
995d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes)
996d71ae5a4SJacob Faibussowitsch {
997bf7f4e0aSPeter Brune   PetscFunctionBegin;
998f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
9996a388c36SPeter Brune   PetscValidPointer(snes, 2);
1000bf7f4e0aSPeter Brune   *snes = linesearch->snes;
10013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1002bf7f4e0aSPeter Brune }
1003bf7f4e0aSPeter Brune 
1004f40b411bSPeter Brune /*@
1005f1c6b773SPeter Brune   SNESLineSearchGetLambda - Gets the last linesearch steplength discovered.
1006f40b411bSPeter Brune 
1007f6dfbefdSBarry Smith   Not Collective
1008f6dfbefdSBarry Smith 
10092fe279fdSBarry Smith   Input Parameter:
10108e557f58SPeter Brune . linesearch - linesearch context
1011f40b411bSPeter Brune 
10122fe279fdSBarry Smith   Output Parameter:
1013f6dfbefdSBarry Smith . lambda - The last steplength computed during `SNESLineSearchApply()`
1014f40b411bSPeter Brune 
101578bcb3b5SPeter Brune   Level: advanced
101678bcb3b5SPeter Brune 
1017f6dfbefdSBarry Smith   Note:
10188e557f58SPeter Brune   This is useful in methods where the solver is ill-scaled and
101978bcb3b5SPeter Brune   requires some adaptive notion of the difference in scale between the
1020f6dfbefdSBarry Smith   solution and the function.  For instance, `SNESQN` may be scaled by the
102178bcb3b5SPeter Brune   line search lambda using the argument -snes_qn_scaling ls.
102278bcb3b5SPeter Brune 
1023f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetLambda()`, `SNESLineSearchGetDamping()`, `SNESLineSearchApply()`
1024f40b411bSPeter Brune @*/
1025d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch linesearch, PetscReal *lambda)
1026d71ae5a4SJacob Faibussowitsch {
10276a388c36SPeter Brune   PetscFunctionBegin;
1028f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1029534a8f05SLisandro Dalcin   PetscValidRealPointer(lambda, 2);
10306a388c36SPeter Brune   *lambda = linesearch->lambda;
10313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10326a388c36SPeter Brune }
10336a388c36SPeter Brune 
1034f40b411bSPeter Brune /*@
1035f6dfbefdSBarry Smith   SNESLineSearchSetLambda - Sets the linesearch steplength
1036f40b411bSPeter Brune 
1037f40b411bSPeter Brune   Input Parameters:
10388e557f58SPeter Brune + linesearch - linesearch context
1039f40b411bSPeter Brune - lambda     - The last steplength.
1040f40b411bSPeter Brune 
104120f4b53cSBarry Smith   Level: advanced
104220f4b53cSBarry Smith 
1043f6dfbefdSBarry Smith   Note:
1044f6dfbefdSBarry Smith   This routine is typically used within implementations of `SNESLineSearchApply()`
1045f6dfbefdSBarry Smith   to set the final steplength.  This routine (and `SNESLineSearchGetLambda()`) were
1046cd7522eaSPeter Brune   added in order to facilitate Quasi-Newton methods that use the previous steplength
1047cd7522eaSPeter Brune   as an inner scaling parameter.
1048cd7522eaSPeter Brune 
1049f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetLambda()`
1050f40b411bSPeter Brune @*/
1051d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
1052d71ae5a4SJacob Faibussowitsch {
10536a388c36SPeter Brune   PetscFunctionBegin;
1054f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
10556a388c36SPeter Brune   linesearch->lambda = lambda;
10563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10576a388c36SPeter Brune }
10586a388c36SPeter Brune 
1059f40b411bSPeter Brune /*@
10603c7d6663SPeter Brune   SNESLineSearchGetTolerances - Gets the tolerances for the linesearch.  These include
106178bcb3b5SPeter Brune   tolerances for the relative and absolute change in the function norm, the change
106278bcb3b5SPeter Brune   in lambda for iterative line searches, the minimum steplength, the maximum steplength,
106378bcb3b5SPeter Brune   and the maximum number of iterations the line search procedure may take.
1064f40b411bSPeter Brune 
1065f6dfbefdSBarry Smith   Not Collective
1066f6dfbefdSBarry Smith 
1067f899ff85SJose E. Roman   Input Parameter:
10688e557f58SPeter Brune . linesearch - linesearch context
1069f40b411bSPeter Brune 
1070f40b411bSPeter Brune   Output Parameters:
1071b13b64c2SBarry Smith + steptol - The minimum steplength
10726cc8e53bSPeter Brune . maxstep - The maximum steplength
1073516fe3c3SPeter Brune . rtol    - The relative tolerance for iterative line searches
1074516fe3c3SPeter Brune . atol    - The absolute tolerance for iterative line searches
1075516fe3c3SPeter Brune . ltol    - The change in lambda tolerance for iterative line searches
1076*e4094ef1SJacob Faibussowitsch - max_its - The maximum number of iterations of the line search
1077f40b411bSPeter Brune 
107878bcb3b5SPeter Brune   Level: intermediate
107978bcb3b5SPeter Brune 
1080f6dfbefdSBarry Smith   Note:
108178bcb3b5SPeter Brune   Different line searches may implement these parameters slightly differently as
10823c7d6663SPeter Brune   the type requires.
1083516fe3c3SPeter Brune 
1084f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetTolerances()`
1085f40b411bSPeter Brune @*/
1086d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch linesearch, PetscReal *steptol, PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
1087d71ae5a4SJacob Faibussowitsch {
10886a388c36SPeter Brune   PetscFunctionBegin;
1089f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1090516fe3c3SPeter Brune   if (steptol) {
1091534a8f05SLisandro Dalcin     PetscValidRealPointer(steptol, 2);
10926a388c36SPeter Brune     *steptol = linesearch->steptol;
1093516fe3c3SPeter Brune   }
1094516fe3c3SPeter Brune   if (maxstep) {
1095534a8f05SLisandro Dalcin     PetscValidRealPointer(maxstep, 3);
1096b13b64c2SBarry Smith     *maxstep = linesearch->maxstep;
1097516fe3c3SPeter Brune   }
1098516fe3c3SPeter Brune   if (rtol) {
1099534a8f05SLisandro Dalcin     PetscValidRealPointer(rtol, 4);
1100516fe3c3SPeter Brune     *rtol = linesearch->rtol;
1101516fe3c3SPeter Brune   }
1102516fe3c3SPeter Brune   if (atol) {
1103534a8f05SLisandro Dalcin     PetscValidRealPointer(atol, 5);
1104516fe3c3SPeter Brune     *atol = linesearch->atol;
1105516fe3c3SPeter Brune   }
1106516fe3c3SPeter Brune   if (ltol) {
1107534a8f05SLisandro Dalcin     PetscValidRealPointer(ltol, 6);
1108516fe3c3SPeter Brune     *ltol = linesearch->ltol;
1109516fe3c3SPeter Brune   }
1110516fe3c3SPeter Brune   if (max_its) {
1111534a8f05SLisandro Dalcin     PetscValidIntPointer(max_its, 7);
1112516fe3c3SPeter Brune     *max_its = linesearch->max_its;
1113516fe3c3SPeter Brune   }
11143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11156a388c36SPeter Brune }
11166a388c36SPeter Brune 
1117f40b411bSPeter Brune /*@
11183c7d6663SPeter Brune   SNESLineSearchSetTolerances -  Gets the tolerances for the linesearch.  These include
111978bcb3b5SPeter Brune   tolerances for the relative and absolute change in the function norm, the change
112078bcb3b5SPeter Brune   in lambda for iterative line searches, the minimum steplength, the maximum steplength,
112178bcb3b5SPeter Brune   and the maximum number of iterations the line search procedure may take.
1122f40b411bSPeter Brune 
1123f6dfbefdSBarry Smith   Collective
1124f6dfbefdSBarry Smith 
1125f40b411bSPeter Brune   Input Parameters:
11268e557f58SPeter Brune + linesearch - linesearch context
1127516fe3c3SPeter Brune . steptol    - The minimum steplength
11286cc8e53bSPeter Brune . maxstep    - The maximum steplength
1129516fe3c3SPeter Brune . rtol       - The relative tolerance for iterative line searches
1130516fe3c3SPeter Brune . atol       - The absolute tolerance for iterative line searches
1131516fe3c3SPeter Brune . ltol       - The change in lambda tolerance for iterative line searches
1132*e4094ef1SJacob Faibussowitsch - max_its    - The maximum number of iterations of the line search
1133f40b411bSPeter Brune 
113420f4b53cSBarry Smith   Level: intermediate
113520f4b53cSBarry Smith 
1136f6dfbefdSBarry Smith   Note:
1137f6dfbefdSBarry Smith   The user may choose to not set any of the tolerances using `PETSC_DEFAULT` in
113878bcb3b5SPeter Brune   place of an argument.
1139f40b411bSPeter Brune 
1140f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetTolerances()`
1141f40b411bSPeter Brune @*/
1142d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch linesearch, PetscReal steptol, PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its)
1143d71ae5a4SJacob Faibussowitsch {
11446a388c36SPeter Brune   PetscFunctionBegin;
1145f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1146d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, steptol, 2);
1147d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, maxstep, 3);
1148d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, rtol, 4);
1149d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, atol, 5);
1150d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, ltol, 6);
1151d3952184SSatish Balay   PetscValidLogicalCollectiveInt(linesearch, max_its, 7);
1152d3952184SSatish Balay 
115313bcc0bdSJacob Faibussowitsch   if (steptol != (PetscReal)PETSC_DEFAULT) {
11545f80ce2aSJacob Faibussowitsch     PetscCheck(steptol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Minimum step length %14.12e must be non-negative", (double)steptol);
11556a388c36SPeter Brune     linesearch->steptol = steptol;
1156d3952184SSatish Balay   }
1157d3952184SSatish Balay 
115813bcc0bdSJacob Faibussowitsch   if (maxstep != (PetscReal)PETSC_DEFAULT) {
11595f80ce2aSJacob Faibussowitsch     PetscCheck(maxstep >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum step length %14.12e must be non-negative", (double)maxstep);
1160516fe3c3SPeter Brune     linesearch->maxstep = maxstep;
1161d3952184SSatish Balay   }
1162d3952184SSatish Balay 
116313bcc0bdSJacob Faibussowitsch   if (rtol != (PetscReal)PETSC_DEFAULT) {
11642061ca32SJunchao 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);
1165516fe3c3SPeter Brune     linesearch->rtol = rtol;
1166d3952184SSatish Balay   }
1167d3952184SSatish Balay 
116813bcc0bdSJacob Faibussowitsch   if (atol != (PetscReal)PETSC_DEFAULT) {
11695f80ce2aSJacob Faibussowitsch     PetscCheck(atol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Absolute tolerance %14.12e must be non-negative", (double)atol);
1170516fe3c3SPeter Brune     linesearch->atol = atol;
1171d3952184SSatish Balay   }
1172d3952184SSatish Balay 
117313bcc0bdSJacob Faibussowitsch   if (ltol != (PetscReal)PETSC_DEFAULT) {
11745f80ce2aSJacob Faibussowitsch     PetscCheck(ltol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Lambda tolerance %14.12e must be non-negative", (double)ltol);
1175516fe3c3SPeter Brune     linesearch->ltol = ltol;
1176d3952184SSatish Balay   }
1177d3952184SSatish Balay 
1178d3952184SSatish Balay   if (max_its != PETSC_DEFAULT) {
117963a3b9bcSJacob Faibussowitsch     PetscCheck(max_its >= 0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of iterations %" PetscInt_FMT " must be non-negative", max_its);
1180516fe3c3SPeter Brune     linesearch->max_its = max_its;
1181d3952184SSatish Balay   }
11823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11836a388c36SPeter Brune }
11846a388c36SPeter Brune 
1185f40b411bSPeter Brune /*@
1186f1c6b773SPeter Brune   SNESLineSearchGetDamping - Gets the line search damping parameter.
1187f40b411bSPeter Brune 
11882fe279fdSBarry Smith   Input Parameter:
11898e557f58SPeter Brune . linesearch - linesearch context
1190f40b411bSPeter Brune 
11912fe279fdSBarry Smith   Output Parameter:
11928e557f58SPeter Brune . damping - The damping parameter
1193f40b411bSPeter Brune 
119478bcb3b5SPeter Brune   Level: advanced
1195f40b411bSPeter Brune 
1196db781477SPatrick Sanan .seealso: `SNESLineSearchGetStepTolerance()`, `SNESQN`
1197f40b411bSPeter Brune @*/
1198d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch linesearch, PetscReal *damping)
1199d71ae5a4SJacob Faibussowitsch {
12006a388c36SPeter Brune   PetscFunctionBegin;
1201f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1202534a8f05SLisandro Dalcin   PetscValidRealPointer(damping, 2);
12036a388c36SPeter Brune   *damping = linesearch->damping;
12043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12056a388c36SPeter Brune }
12066a388c36SPeter Brune 
1207f40b411bSPeter Brune /*@
1208fd292e60Sprj-   SNESLineSearchSetDamping - Sets the line search damping parameter.
1209f40b411bSPeter Brune 
1210f40b411bSPeter Brune   Input Parameters:
121103fd4120SBarry Smith + linesearch - linesearch context
121203fd4120SBarry Smith - damping    - The damping parameter
1213f40b411bSPeter Brune 
1214f6dfbefdSBarry Smith   Options Database Key:
1215f6dfbefdSBarry Smith . -snes_linesearch_damping <damping> - the damping value
1216f6dfbefdSBarry Smith 
1217f40b411bSPeter Brune   Level: intermediate
1218f40b411bSPeter Brune 
1219f6dfbefdSBarry Smith   Note:
1220f6dfbefdSBarry Smith   The `SNESLINESEARCHNONE` line search merely takes the update step scaled by the damping parameter.
1221cd7522eaSPeter Brune   The use of the damping parameter in the l2 and cp line searches is much more subtle;
122278bcb3b5SPeter Brune   it is used as a starting point in calculating the secant step. However, the eventual
1223cd7522eaSPeter Brune   step may be of greater length than the damping parameter.  In the bt line search it is
1224cd7522eaSPeter Brune   used as the maximum possible step length, as the bt line search only backtracks.
1225cd7522eaSPeter Brune 
1226f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetDamping()`
1227f40b411bSPeter Brune @*/
1228d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch linesearch, PetscReal damping)
1229d71ae5a4SJacob Faibussowitsch {
12306a388c36SPeter Brune   PetscFunctionBegin;
1231f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
12326a388c36SPeter Brune   linesearch->damping = damping;
12333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12346a388c36SPeter Brune }
12356a388c36SPeter Brune 
123659405d5eSPeter Brune /*@
123759405d5eSPeter Brune   SNESLineSearchGetOrder - Gets the line search approximation order.
123859405d5eSPeter Brune 
1239f6dfbefdSBarry Smith   Input Parameter:
124078bcb3b5SPeter Brune . linesearch - linesearch context
124159405d5eSPeter Brune 
1242f6dfbefdSBarry Smith   Output Parameter:
12438e557f58SPeter Brune . order - The order
124459405d5eSPeter Brune 
124578bcb3b5SPeter Brune   Possible Values for order:
1246f6dfbefdSBarry Smith +  1 or `SNES_LINESEARCH_ORDER_LINEAR` - linear order
1247f6dfbefdSBarry Smith .  2 or `SNES_LINESEARCH_ORDER_QUADRATIC` - quadratic order
1248f6dfbefdSBarry Smith -  3 or `SNES_LINESEARCH_ORDER_CUBIC` - cubic order
124978bcb3b5SPeter Brune 
125059405d5eSPeter Brune   Level: intermediate
125159405d5eSPeter Brune 
1252f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetOrder()`
125359405d5eSPeter Brune @*/
1254d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch linesearch, PetscInt *order)
1255d71ae5a4SJacob Faibussowitsch {
125659405d5eSPeter Brune   PetscFunctionBegin;
125759405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1258534a8f05SLisandro Dalcin   PetscValidIntPointer(order, 2);
125959405d5eSPeter Brune   *order = linesearch->order;
12603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
126159405d5eSPeter Brune }
126259405d5eSPeter Brune 
126359405d5eSPeter Brune /*@
12641f8196a2SJed Brown   SNESLineSearchSetOrder - Sets the maximum order of the polynomial fit used in the line search
126559405d5eSPeter Brune 
126659405d5eSPeter Brune   Input Parameters:
1267a2b725a8SWilliam Gropp + linesearch - linesearch context
1268a2b725a8SWilliam Gropp - order      - The damping parameter
126959405d5eSPeter Brune 
127059405d5eSPeter Brune   Level: intermediate
127159405d5eSPeter Brune 
127278bcb3b5SPeter Brune   Possible Values for order:
1273f6dfbefdSBarry Smith +  1 or `SNES_LINESEARCH_ORDER_LINEAR` - linear order
1274f6dfbefdSBarry Smith .  2 or `SNES_LINESEARCH_ORDER_QUADRATIC` - quadratic order
1275f6dfbefdSBarry Smith -  3 or `SNES_LINESEARCH_ORDER_CUBIC` - cubic order
127678bcb3b5SPeter Brune 
127759405d5eSPeter Brune   Notes:
127859405d5eSPeter Brune   Variable orders are supported by the following line searches:
127978bcb3b5SPeter Brune +  bt - cubic and quadratic
128078bcb3b5SPeter Brune -  cp - linear and quadratic
128159405d5eSPeter Brune 
1282f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetOrder()`, `SNESLineSearchSetDamping()`
128359405d5eSPeter Brune @*/
1284d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch linesearch, PetscInt order)
1285d71ae5a4SJacob Faibussowitsch {
128659405d5eSPeter Brune   PetscFunctionBegin;
128759405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
128859405d5eSPeter Brune   linesearch->order = order;
12893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
129059405d5eSPeter Brune }
129159405d5eSPeter Brune 
1292f40b411bSPeter Brune /*@
1293f1c6b773SPeter Brune   SNESLineSearchGetNorms - Gets the norms for for X, Y, and F.
1294f40b411bSPeter Brune 
1295f6dfbefdSBarry Smith   Not Collective
1296f6dfbefdSBarry Smith 
1297f899ff85SJose E. Roman   Input Parameter:
129878bcb3b5SPeter Brune . linesearch - linesearch context
1299f40b411bSPeter Brune 
1300f40b411bSPeter Brune   Output Parameters:
1301f40b411bSPeter Brune + xnorm - The norm of the current solution
1302a68bbae5SBarry Smith . fnorm - The norm of the current function, this is the `norm(function(X))` where `X` is the current solution.
1303f40b411bSPeter Brune - ynorm - The norm of the current update
1304f40b411bSPeter Brune 
130578bcb3b5SPeter Brune   Level: developer
1306f40b411bSPeter Brune 
1307a68bbae5SBarry Smith   Note:
1308a68bbae5SBarry Smith   Some values may not be up to date at particular points in the code.
1309a68bbae5SBarry Smith 
1310a68bbae5SBarry Smith   This, in combination with `SNESLineSearchSetNorms()`, allow the line search and the `SNESSolve_XXX()` to share
1311a68bbae5SBarry Smith   computed values.
1312a68bbae5SBarry Smith 
1313f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetNorms()` `SNESLineSearchGetVecs()`
1314f40b411bSPeter Brune @*/
1315d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal *xnorm, PetscReal *fnorm, PetscReal *ynorm)
1316d71ae5a4SJacob Faibussowitsch {
1317bf7f4e0aSPeter Brune   PetscFunctionBegin;
1318f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1319f5af7f23SKarl Rupp   if (xnorm) *xnorm = linesearch->xnorm;
1320f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
1321f5af7f23SKarl Rupp   if (ynorm) *ynorm = linesearch->ynorm;
13223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1323bf7f4e0aSPeter Brune }
1324bf7f4e0aSPeter Brune 
1325f40b411bSPeter Brune /*@
1326f1c6b773SPeter Brune   SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F.
1327f40b411bSPeter Brune 
1328c3339decSBarry Smith   Collective
1329f6dfbefdSBarry Smith 
1330f40b411bSPeter Brune   Input Parameters:
133178bcb3b5SPeter Brune + linesearch - linesearch context
1332f40b411bSPeter Brune . xnorm      - The norm of the current solution
1333a68bbae5SBarry Smith . fnorm      - The norm of the current function, this is the `norm(function(X))` where `X` is the current solution
1334f40b411bSPeter Brune - ynorm      - The norm of the current update
1335f40b411bSPeter Brune 
1336f6dfbefdSBarry Smith   Level: developer
1337f40b411bSPeter Brune 
1338f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()`
1339f40b411bSPeter Brune @*/
1340d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
1341d71ae5a4SJacob Faibussowitsch {
13426a388c36SPeter Brune   PetscFunctionBegin;
1343f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
13446a388c36SPeter Brune   linesearch->xnorm = xnorm;
13456a388c36SPeter Brune   linesearch->fnorm = fnorm;
13466a388c36SPeter Brune   linesearch->ynorm = ynorm;
13473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13486a388c36SPeter Brune }
13496a388c36SPeter Brune 
1350f40b411bSPeter Brune /*@
1351f1c6b773SPeter Brune   SNESLineSearchComputeNorms - Computes the norms of X, F, and Y.
1352f40b411bSPeter Brune 
13532fe279fdSBarry Smith   Input Parameter:
135478bcb3b5SPeter Brune . linesearch - linesearch context
1355f40b411bSPeter Brune 
135620f4b53cSBarry Smith   Options Database Key:
13578e557f58SPeter Brune . -snes_linesearch_norms - turn norm computation on or off
1358f40b411bSPeter Brune 
1359f40b411bSPeter Brune   Level: intermediate
1360f40b411bSPeter Brune 
1361f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetNorms`, `SNESLineSearchSetNorms()`, `SNESLineSearchSetComputeNorms()`
1362f40b411bSPeter Brune @*/
1363d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
1364d71ae5a4SJacob Faibussowitsch {
13659bd66eb0SPeter Brune   SNES snes;
1366bf388a1fSBarry Smith 
13676a388c36SPeter Brune   PetscFunctionBegin;
13686a388c36SPeter Brune   if (linesearch->norms) {
13699bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
13709566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
13719566063dSJacob Faibussowitsch       PetscCall(VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
13729566063dSJacob Faibussowitsch       PetscCall(VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm));
13739566063dSJacob Faibussowitsch       PetscCall((*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm));
13749bd66eb0SPeter Brune     } else {
13759566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(linesearch->vec_func, NORM_2, &linesearch->fnorm));
13769566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
13779566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm));
13789566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(linesearch->vec_func, NORM_2, &linesearch->fnorm));
13799566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
13809566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(linesearch->vec_update, NORM_2, &linesearch->ynorm));
13816a388c36SPeter Brune     }
13829bd66eb0SPeter Brune   }
13833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13846a388c36SPeter Brune }
13856a388c36SPeter Brune 
13866f263ca3SPeter Brune /*@
13876f263ca3SPeter Brune   SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search.
13886f263ca3SPeter Brune 
13896f263ca3SPeter Brune   Input Parameters:
139078bcb3b5SPeter Brune + linesearch - linesearch context
139178bcb3b5SPeter Brune - flg        - indicates whether or not to compute norms
13926f263ca3SPeter Brune 
139320f4b53cSBarry Smith   Options Database Key:
13940b00b554SBarry Smith . -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic (none) linesearch
13956f263ca3SPeter Brune 
139620f4b53cSBarry Smith   Level: intermediate
139720f4b53cSBarry Smith 
1398f6dfbefdSBarry Smith   Note:
1399f6dfbefdSBarry 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.
14006f263ca3SPeter Brune 
1401f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetNorms()`, `SNESLineSearchComputeNorms()`, `SNESLINESEARCHBASIC`
14026f263ca3SPeter Brune @*/
1403d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg)
1404d71ae5a4SJacob Faibussowitsch {
14056f263ca3SPeter Brune   PetscFunctionBegin;
14066f263ca3SPeter Brune   linesearch->norms = flg;
14073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14086f263ca3SPeter Brune }
14096f263ca3SPeter Brune 
1410f40b411bSPeter Brune /*@
1411f6dfbefdSBarry Smith   SNESLineSearchGetVecs - Gets the vectors from the `SNESLineSearch` context
1412f6dfbefdSBarry Smith 
1413f6dfbefdSBarry Smith   Not Collective on but the vectors are parallel
1414f40b411bSPeter Brune 
1415f899ff85SJose E. Roman   Input Parameter:
141678bcb3b5SPeter Brune . linesearch - linesearch context
1417f40b411bSPeter Brune 
1418f40b411bSPeter Brune   Output Parameters:
14196232e825SPeter Brune + X - Solution vector
14206232e825SPeter Brune . F - Function vector
14216232e825SPeter Brune . Y - Search direction vector
14226232e825SPeter Brune . W - Solution work vector
14236232e825SPeter Brune - G - Function work vector
14246232e825SPeter Brune 
142520f4b53cSBarry Smith   Level: advanced
142620f4b53cSBarry Smith 
14277bba9028SPeter Brune   Notes:
142820f4b53cSBarry Smith   At the beginning of a line search application, `X` should contain a
142920f4b53cSBarry Smith   solution and the vector `F` the function computed at `X`.  At the end of the
143020f4b53cSBarry Smith   line search application, `X` should contain the new solution, and `F` the
14316232e825SPeter Brune   function evaluated at the new solution.
1432f40b411bSPeter Brune 
1433f6dfbefdSBarry Smith   These vectors are owned by the `SNESLineSearch` and should not be destroyed by the caller
14342a7a6963SBarry Smith 
1435f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()`
1436f40b411bSPeter Brune @*/
1437d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch, Vec *X, Vec *F, Vec *Y, Vec *W, Vec *G)
1438d71ae5a4SJacob Faibussowitsch {
14396a388c36SPeter Brune   PetscFunctionBegin;
1440f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
14416a388c36SPeter Brune   if (X) {
14426a388c36SPeter Brune     PetscValidPointer(X, 2);
14436a388c36SPeter Brune     *X = linesearch->vec_sol;
14446a388c36SPeter Brune   }
14456a388c36SPeter Brune   if (F) {
14466a388c36SPeter Brune     PetscValidPointer(F, 3);
14476a388c36SPeter Brune     *F = linesearch->vec_func;
14486a388c36SPeter Brune   }
14496a388c36SPeter Brune   if (Y) {
14506a388c36SPeter Brune     PetscValidPointer(Y, 4);
14516a388c36SPeter Brune     *Y = linesearch->vec_update;
14526a388c36SPeter Brune   }
14536a388c36SPeter Brune   if (W) {
14546a388c36SPeter Brune     PetscValidPointer(W, 5);
14556a388c36SPeter Brune     *W = linesearch->vec_sol_new;
14566a388c36SPeter Brune   }
14576a388c36SPeter Brune   if (G) {
14586a388c36SPeter Brune     PetscValidPointer(G, 6);
14596a388c36SPeter Brune     *G = linesearch->vec_func_new;
14606a388c36SPeter Brune   }
14613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14626a388c36SPeter Brune }
14636a388c36SPeter Brune 
1464f40b411bSPeter Brune /*@
1465f6dfbefdSBarry Smith   SNESLineSearchSetVecs - Sets the vectors on the `SNESLineSearch` context
1466f6dfbefdSBarry Smith 
1467c3339decSBarry Smith   Logically Collective
1468f40b411bSPeter Brune 
1469f40b411bSPeter Brune   Input Parameters:
147078bcb3b5SPeter Brune + linesearch - linesearch context
14716232e825SPeter Brune . X          - Solution vector
14726232e825SPeter Brune . F          - Function vector
14736232e825SPeter Brune . Y          - Search direction vector
14746232e825SPeter Brune . W          - Solution work vector
14756232e825SPeter Brune - G          - Function work vector
1476f40b411bSPeter Brune 
147778bcb3b5SPeter Brune   Level: advanced
1478f40b411bSPeter Brune 
1479f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetNorms()`, `SNESLineSearchGetVecs()`
1480f40b411bSPeter Brune @*/
1481d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch, Vec X, Vec F, Vec Y, Vec W, Vec G)
1482d71ae5a4SJacob Faibussowitsch {
14836a388c36SPeter Brune   PetscFunctionBegin;
1484f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
14856a388c36SPeter Brune   if (X) {
14866a388c36SPeter Brune     PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
14876a388c36SPeter Brune     linesearch->vec_sol = X;
14886a388c36SPeter Brune   }
14896a388c36SPeter Brune   if (F) {
14906a388c36SPeter Brune     PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
14916a388c36SPeter Brune     linesearch->vec_func = F;
14926a388c36SPeter Brune   }
14936a388c36SPeter Brune   if (Y) {
14946a388c36SPeter Brune     PetscValidHeaderSpecific(Y, VEC_CLASSID, 4);
14956a388c36SPeter Brune     linesearch->vec_update = Y;
14966a388c36SPeter Brune   }
14976a388c36SPeter Brune   if (W) {
14986a388c36SPeter Brune     PetscValidHeaderSpecific(W, VEC_CLASSID, 5);
14996a388c36SPeter Brune     linesearch->vec_sol_new = W;
15006a388c36SPeter Brune   }
15016a388c36SPeter Brune   if (G) {
15026a388c36SPeter Brune     PetscValidHeaderSpecific(G, VEC_CLASSID, 6);
15036a388c36SPeter Brune     linesearch->vec_func_new = G;
15046a388c36SPeter Brune   }
15053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15066a388c36SPeter Brune }
15076a388c36SPeter Brune 
1508e7058c64SPeter Brune /*@C
1509f1c6b773SPeter Brune   SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1510f6dfbefdSBarry Smith   `SNESLineSearch` options in the database.
1511e7058c64SPeter Brune 
1512c3339decSBarry Smith   Logically Collective
1513e7058c64SPeter Brune 
1514e7058c64SPeter Brune   Input Parameters:
1515f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context
1516e7058c64SPeter Brune - prefix     - the prefix to prepend to all option names
1517e7058c64SPeter Brune 
151820f4b53cSBarry Smith   Level: advanced
151920f4b53cSBarry Smith 
1520f6dfbefdSBarry Smith   Note:
1521e7058c64SPeter Brune   A hyphen (-) must NOT be given at the beginning of the prefix name.
1522e7058c64SPeter Brune   The first character of all runtime options is AUTOMATICALLY the hyphen.
1523e7058c64SPeter Brune 
1524f6dfbefdSBarry Smith .seealso: `SNESLineSearch()`, `SNESLineSearchSetFromOptions()`, `SNESGetOptionsPrefix()`
1525e7058c64SPeter Brune @*/
1526d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch, const char prefix[])
1527d71ae5a4SJacob Faibussowitsch {
1528e7058c64SPeter Brune   PetscFunctionBegin;
1529f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
15309566063dSJacob Faibussowitsch   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)linesearch, prefix));
15313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1532e7058c64SPeter Brune }
1533e7058c64SPeter Brune 
1534e7058c64SPeter Brune /*@C
1535f6dfbefdSBarry Smith   SNESLineSearchGetOptionsPrefix - Gets the prefix used for searching for all
1536f1c6b773SPeter Brune   SNESLineSearch options in the database.
1537e7058c64SPeter Brune 
1538e7058c64SPeter Brune   Not Collective
1539e7058c64SPeter Brune 
1540e7058c64SPeter Brune   Input Parameter:
1541f6dfbefdSBarry Smith . linesearch - the `SNESLineSearch` context
1542e7058c64SPeter Brune 
1543e7058c64SPeter Brune   Output Parameter:
1544e7058c64SPeter Brune . prefix - pointer to the prefix string used
1545e7058c64SPeter Brune 
1546e7058c64SPeter Brune   Level: advanced
1547e7058c64SPeter Brune 
1548*e4094ef1SJacob Faibussowitsch   Fortran Notes:
154920f4b53cSBarry Smith   The user should pass in a string 'prefix' of
155020f4b53cSBarry Smith   sufficient length to hold the prefix.
155120f4b53cSBarry Smith 
1552f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESAppendOptionsPrefix()`
1553e7058c64SPeter Brune @*/
1554d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch, const char *prefix[])
1555d71ae5a4SJacob Faibussowitsch {
1556e7058c64SPeter Brune   PetscFunctionBegin;
1557f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
15589566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)linesearch, prefix));
15593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1560e7058c64SPeter Brune }
1561bf7f4e0aSPeter Brune 
15628d359177SBarry Smith /*@C
1563f6dfbefdSBarry Smith   SNESLineSearchSetWorkVecs - Sets work vectors for the line search.
1564f40b411bSPeter Brune 
1565d8d19677SJose E. Roman   Input Parameters:
1566f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context
1567f40b411bSPeter Brune - nwork      - the number of work vectors
1568f40b411bSPeter Brune 
1569f40b411bSPeter Brune   Level: developer
1570f40b411bSPeter Brune 
1571f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESSetWorkVecs()`
1572f40b411bSPeter Brune @*/
1573d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork)
1574d71ae5a4SJacob Faibussowitsch {
1575bf7f4e0aSPeter Brune   PetscFunctionBegin;
15760fdf79fbSJacob Faibussowitsch   PetscCheck(linesearch->vec_sol, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
15779566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work));
15783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1579bf7f4e0aSPeter Brune }
1580bf7f4e0aSPeter Brune 
1581f40b411bSPeter Brune /*@
1582422a814eSBarry Smith   SNESLineSearchGetReason - Gets the success/failure status of the last line search application
1583f40b411bSPeter Brune 
15842fe279fdSBarry Smith   Input Parameter:
158578bcb3b5SPeter Brune . linesearch - linesearch context
1586f40b411bSPeter Brune 
15872fe279fdSBarry Smith   Output Parameter:
1588422a814eSBarry Smith . result - The success or failure status
1589f40b411bSPeter Brune 
159020f4b53cSBarry Smith   Level: developer
159120f4b53cSBarry Smith 
1592f6dfbefdSBarry Smith   Note:
1593f6dfbefdSBarry Smith   This is typically called after `SNESLineSearchApply()` in order to determine if the line-search failed
1594cd7522eaSPeter Brune   (and set the SNES convergence accordingly).
1595cd7522eaSPeter Brune 
1596f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetReason()`, `SNESLineSearchReason`
1597f40b411bSPeter Brune @*/
1598d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result)
1599d71ae5a4SJacob Faibussowitsch {
1600bf7f4e0aSPeter Brune   PetscFunctionBegin;
1601f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1602422a814eSBarry Smith   PetscValidPointer(result, 2);
1603422a814eSBarry Smith   *result = linesearch->result;
16043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1605bf7f4e0aSPeter Brune }
1606bf7f4e0aSPeter Brune 
1607f40b411bSPeter Brune /*@
1608422a814eSBarry Smith   SNESLineSearchSetReason - Sets the success/failure status of the last line search application
1609f40b411bSPeter Brune 
1610f40b411bSPeter Brune   Input Parameters:
161178bcb3b5SPeter Brune + linesearch - linesearch context
1612422a814eSBarry Smith - result     - The success or failure status
1613f40b411bSPeter Brune 
161420f4b53cSBarry Smith   Level: developer
161520f4b53cSBarry Smith 
1616f6dfbefdSBarry Smith   Note:
1617f6dfbefdSBarry Smith   This is typically called in a `SNESLineSearchApply()` or a `SNESLINESEARCHSHELL` implementation to set
1618cd7522eaSPeter Brune   the success or failure of the line search method.
1619cd7522eaSPeter Brune 
1620f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchReason`, `SNESLineSearchGetSResult()`
1621f40b411bSPeter Brune @*/
1622d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result)
1623d71ae5a4SJacob Faibussowitsch {
16246a388c36SPeter Brune   PetscFunctionBegin;
16255fd66863SKarl Rupp   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1626422a814eSBarry Smith   linesearch->result = result;
16273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16286a388c36SPeter Brune }
16296a388c36SPeter Brune 
16309bd66eb0SPeter Brune /*@C
1631f1c6b773SPeter Brune   SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.
16329bd66eb0SPeter Brune 
1633c3339decSBarry Smith   Logically Collective
1634f6dfbefdSBarry Smith 
16359bd66eb0SPeter Brune   Input Parameters:
1636f6dfbefdSBarry Smith +  snes - nonlinear context obtained from `SNESCreate()`
16379bd66eb0SPeter Brune . projectfunc - function for projecting the function to the bounds
16389bd66eb0SPeter Brune - normfunc    - function for computing the norm of an active set
16399bd66eb0SPeter Brune 
164020f4b53cSBarry Smith   Calling sequence of `projectfunc`:
16419bd66eb0SPeter Brune .vb
164220f4b53cSBarry Smith    PetscErrorCode projectfunc(SNES snes, Vec X)
16439bd66eb0SPeter Brune .ve
16449bd66eb0SPeter Brune +   snes - nonlinear context
164520f4b53cSBarry Smith -   X - current solution, store the projected solution here
16469bd66eb0SPeter Brune 
164720f4b53cSBarry Smith   Calling sequence of `normfunc`:
16489bd66eb0SPeter Brune .vb
164920f4b53cSBarry Smith    PetscErrorCode normfunc(SNES snes, Vec X, Vec F, PetscScalar *fnorm)
16509bd66eb0SPeter Brune .ve
16519bd66eb0SPeter Brune +   snes - nonlinear context
16529bd66eb0SPeter Brune .   X - current solution
165320f4b53cSBarry Smith .   F - current residual
1654*e4094ef1SJacob Faibussowitsch - linesearch - VI-specific norm of the function
16559bd66eb0SPeter Brune 
1656f6dfbefdSBarry Smith   Level: advanced
16579bd66eb0SPeter Brune 
165820f4b53cSBarry Smith   Note:
165920f4b53cSBarry Smith   The VI solvers require projection of the solution to the feasible set.  `projectfunc` should implement this.
166020f4b53cSBarry Smith 
166120f4b53cSBarry Smith   The VI solvers require special evaluation of the function norm such that the norm is only calculated
166220f4b53cSBarry Smith   on the inactive set.  This should be implemented by `normfunc`.
166320f4b53cSBarry Smith 
1664f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchGetVIFunctions()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()`
16659bd66eb0SPeter Brune @*/
1666d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc)
1667d71ae5a4SJacob Faibussowitsch {
16689bd66eb0SPeter Brune   PetscFunctionBegin;
1669f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
16709bd66eb0SPeter Brune   if (projectfunc) linesearch->ops->viproject = projectfunc;
16719bd66eb0SPeter Brune   if (normfunc) linesearch->ops->vinorm = normfunc;
16723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16739bd66eb0SPeter Brune }
16749bd66eb0SPeter Brune 
16759bd66eb0SPeter Brune /*@C
1676f1c6b773SPeter Brune   SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.
16779bd66eb0SPeter Brune 
1678f6dfbefdSBarry Smith   Not Collective
1679f6dfbefdSBarry Smith 
1680f899ff85SJose E. Roman   Input Parameter:
1681f6dfbefdSBarry Smith . linesearch - the line search context, obtain with `SNESGetLineSearch()`
16829bd66eb0SPeter Brune 
16839bd66eb0SPeter Brune   Output Parameters:
16849bd66eb0SPeter Brune + projectfunc - function for projecting the function to the bounds
16859bd66eb0SPeter Brune - normfunc    - function for computing the norm of an active set
16869bd66eb0SPeter Brune 
1687f6dfbefdSBarry Smith   Level: advanced
16889bd66eb0SPeter Brune 
1689f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchSetVIFunctions()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()`
16909bd66eb0SPeter Brune @*/
1691d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc)
1692d71ae5a4SJacob Faibussowitsch {
16939bd66eb0SPeter Brune   PetscFunctionBegin;
16949bd66eb0SPeter Brune   if (projectfunc) *projectfunc = linesearch->ops->viproject;
16959bd66eb0SPeter Brune   if (normfunc) *normfunc = linesearch->ops->vinorm;
16963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16979bd66eb0SPeter Brune }
16989bd66eb0SPeter Brune 
1699bf7f4e0aSPeter Brune /*@C
1700f6dfbefdSBarry Smith   SNESLineSearchRegister - register a line search method
1701bf7f4e0aSPeter Brune 
1702bf7f4e0aSPeter Brune   Level: advanced
1703f6dfbefdSBarry Smith 
1704f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchSetType()`
1705bf7f4e0aSPeter Brune @*/
1706d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchRegister(const char sname[], PetscErrorCode (*function)(SNESLineSearch))
1707d71ae5a4SJacob Faibussowitsch {
1708bf7f4e0aSPeter Brune   PetscFunctionBegin;
17099566063dSJacob Faibussowitsch   PetscCall(SNESInitializePackage());
17109566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&SNESLineSearchList, sname, function));
17113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1712bf7f4e0aSPeter Brune }
1713