xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 /*@
10dcf2fd19SBarry Smith    SNESLineSearchMonitorCancel - Clears all the monitor functions for a SNESLineSearch object.
11dcf2fd19SBarry Smith 
12dcf2fd19SBarry Smith    Logically Collective on SNESLineSearch
13dcf2fd19SBarry Smith 
14dcf2fd19SBarry Smith    Input Parameters:
15dcf2fd19SBarry 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
19dcf2fd19SBarry Smith     into a code by calls to SNESLineSearchMonitorSet(), but does not cancel those
20dcf2fd19SBarry Smith     set via the options database
21dcf2fd19SBarry Smith 
22dcf2fd19SBarry Smith    Notes:
23dcf2fd19SBarry Smith    There is no way to clear one specific monitor from a SNESLineSearch object.
24dcf2fd19SBarry Smith 
25dcf2fd19SBarry Smith    This does not clear the monitor set with SNESLineSearchSetDefaultMonitor() use SNESLineSearchSetDefaultMonitor(ls,NULL) to cancel
26dcf2fd19SBarry Smith    that one.
27dcf2fd19SBarry Smith 
28dcf2fd19SBarry Smith    Level: intermediate
29dcf2fd19SBarry Smith 
30db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchMonitorDefault()`, `SNESLineSearchMonitorSet()`
31dcf2fd19SBarry Smith @*/
329371c9d4SSatish Balay PetscErrorCode SNESLineSearchMonitorCancel(SNESLineSearch ls) {
33dcf2fd19SBarry Smith   PetscInt i;
34dcf2fd19SBarry Smith 
35dcf2fd19SBarry Smith   PetscFunctionBegin;
36dcf2fd19SBarry Smith   PetscValidHeaderSpecific(ls, SNESLINESEARCH_CLASSID, 1);
37dcf2fd19SBarry Smith   for (i = 0; i < ls->numbermonitors; i++) {
38*48a46eb9SPierre Jolivet     if (ls->monitordestroy[i]) PetscCall((*ls->monitordestroy[i])(&ls->monitorcontext[i]));
39dcf2fd19SBarry Smith   }
40dcf2fd19SBarry Smith   ls->numbermonitors = 0;
41dcf2fd19SBarry Smith   PetscFunctionReturn(0);
42dcf2fd19SBarry Smith }
43dcf2fd19SBarry Smith 
44dcf2fd19SBarry Smith /*@
45dcf2fd19SBarry Smith    SNESLineSearchMonitor - runs the user provided monitor routines, if they exist
46dcf2fd19SBarry Smith 
47dcf2fd19SBarry Smith    Collective on SNES
48dcf2fd19SBarry Smith 
49dcf2fd19SBarry Smith    Input Parameters:
50dcf2fd19SBarry Smith .  ls - the linesearch object
51dcf2fd19SBarry Smith 
52dcf2fd19SBarry Smith    Notes:
53dcf2fd19SBarry Smith    This routine is called by the SNES implementations.
54dcf2fd19SBarry Smith    It does not typically need to be called by the user.
55dcf2fd19SBarry Smith 
56dcf2fd19SBarry Smith    Level: developer
57dcf2fd19SBarry Smith 
58db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchMonitorSet()`
59dcf2fd19SBarry Smith @*/
609371c9d4SSatish Balay PetscErrorCode SNESLineSearchMonitor(SNESLineSearch ls) {
61dcf2fd19SBarry Smith   PetscInt i, n = ls->numbermonitors;
62dcf2fd19SBarry Smith 
63dcf2fd19SBarry Smith   PetscFunctionBegin;
64*48a46eb9SPierre Jolivet   for (i = 0; i < n; i++) PetscCall((*ls->monitorftns[i])(ls, ls->monitorcontext[i]));
65dcf2fd19SBarry Smith   PetscFunctionReturn(0);
66dcf2fd19SBarry Smith }
67dcf2fd19SBarry Smith 
68dcf2fd19SBarry Smith /*@C
69dcf2fd19SBarry Smith    SNESLineSearchMonitorSet - Sets an ADDITIONAL function that is to be used at every
70dcf2fd19SBarry Smith    iteration of the nonlinear solver to display the iteration's
71dcf2fd19SBarry Smith    progress.
72dcf2fd19SBarry Smith 
73dcf2fd19SBarry Smith    Logically Collective on SNESLineSearch
74dcf2fd19SBarry Smith 
75dcf2fd19SBarry Smith    Input Parameters:
76dcf2fd19SBarry Smith +  ls - the SNESLineSearch context
77dcf2fd19SBarry Smith .  f - the monitor function
78dcf2fd19SBarry Smith .  mctx - [optional] user-defined context for private data for the
79dcf2fd19SBarry Smith           monitor routine (use NULL if no context is desired)
80dcf2fd19SBarry Smith -  monitordestroy - [optional] routine that frees monitor context
81dcf2fd19SBarry Smith           (may be NULL)
82dcf2fd19SBarry Smith 
83dcf2fd19SBarry Smith    Notes:
84dcf2fd19SBarry Smith    Several different monitoring routines may be set by calling
85dcf2fd19SBarry Smith    SNESLineSearchMonitorSet() multiple times; all will be called in the
86dcf2fd19SBarry Smith    order in which they were set.
87dcf2fd19SBarry Smith 
8895452b02SPatrick Sanan    Fortran Notes:
8995452b02SPatrick Sanan     Only a single monitor function can be set for each SNESLineSearch object
90dcf2fd19SBarry Smith 
91dcf2fd19SBarry Smith    Level: intermediate
92dcf2fd19SBarry Smith 
93db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchMonitorDefault()`, `SNESLineSearchMonitorCancel()`
94dcf2fd19SBarry Smith @*/
959371c9d4SSatish Balay PetscErrorCode SNESLineSearchMonitorSet(SNESLineSearch ls, PetscErrorCode (*f)(SNESLineSearch, void *), void *mctx, PetscErrorCode (*monitordestroy)(void **)) {
9678064530SBarry Smith   PetscInt  i;
9778064530SBarry Smith   PetscBool identical;
9878064530SBarry Smith 
99dcf2fd19SBarry Smith   PetscFunctionBegin;
100dcf2fd19SBarry Smith   PetscValidHeaderSpecific(ls, SNESLINESEARCH_CLASSID, 1);
10178064530SBarry Smith   for (i = 0; i < ls->numbermonitors; i++) {
1029566063dSJacob Faibussowitsch     PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))f, mctx, monitordestroy, (PetscErrorCode(*)(void))ls->monitorftns[i], ls->monitorcontext[i], ls->monitordestroy[i], &identical));
10378064530SBarry Smith     if (identical) PetscFunctionReturn(0);
10478064530SBarry Smith   }
1055f80ce2aSJacob Faibussowitsch   PetscCheck(ls->numbermonitors < MAXSNESLSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many monitors set");
106dcf2fd19SBarry Smith   ls->monitorftns[ls->numbermonitors]      = f;
107dcf2fd19SBarry Smith   ls->monitordestroy[ls->numbermonitors]   = monitordestroy;
108dcf2fd19SBarry Smith   ls->monitorcontext[ls->numbermonitors++] = (void *)mctx;
109dcf2fd19SBarry Smith   PetscFunctionReturn(0);
110dcf2fd19SBarry Smith }
111dcf2fd19SBarry Smith 
112dcf2fd19SBarry Smith /*@C
113dcf2fd19SBarry Smith    SNESLineSearchMonitorSolutionUpdate - Monitors each update a new function value the linesearch tries
114dcf2fd19SBarry Smith 
115dcf2fd19SBarry Smith    Collective on SNESLineSearch
116dcf2fd19SBarry Smith 
117dcf2fd19SBarry Smith    Input Parameters:
118dcf2fd19SBarry Smith +  ls - the SNES linesearch object
119d12e167eSBarry Smith -  vf - the context for the monitor, in this case it is an ASCII PetscViewer and format
120dcf2fd19SBarry Smith 
121dcf2fd19SBarry Smith    Level: intermediate
122dcf2fd19SBarry Smith 
123db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESMonitorSet()`, `SNESMonitorSolution()`
124dcf2fd19SBarry Smith @*/
1259371c9d4SSatish Balay PetscErrorCode SNESLineSearchMonitorSolutionUpdate(SNESLineSearch ls, PetscViewerAndFormat *vf) {
126d12e167eSBarry Smith   PetscViewer viewer = vf->viewer;
127dcf2fd19SBarry Smith   Vec         Y, W, G;
128dcf2fd19SBarry Smith 
129dcf2fd19SBarry Smith   PetscFunctionBegin;
1309566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetVecs(ls, NULL, NULL, &Y, &W, &G));
1319566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(viewer, vf->format));
1329566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted update to solution \n"));
1339566063dSJacob Faibussowitsch   PetscCall(VecView(Y, viewer));
1349566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted new solution \n"));
1359566063dSJacob Faibussowitsch   PetscCall(VecView(W, viewer));
1369566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted updated function value\n"));
1379566063dSJacob Faibussowitsch   PetscCall(VecView(G, viewer));
1389566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(viewer));
139dcf2fd19SBarry Smith   PetscFunctionReturn(0);
140dcf2fd19SBarry Smith }
141dcf2fd19SBarry Smith 
142f40b411bSPeter Brune /*@
143cd7522eaSPeter Brune    SNESLineSearchCreate - Creates the line search context.
144f40b411bSPeter Brune 
145cd7522eaSPeter Brune    Logically Collective on Comm
146f40b411bSPeter Brune 
147f40b411bSPeter Brune    Input Parameters:
148cd7522eaSPeter Brune .  comm - MPI communicator for the line search (typically from the associated SNES context).
149f40b411bSPeter Brune 
150f40b411bSPeter Brune    Output Parameters:
1518e557f58SPeter Brune .  outlinesearch - the new linesearch context
152f40b411bSPeter Brune 
153162e0bf5SPeter Brune    Level: developer
154162e0bf5SPeter Brune 
155162e0bf5SPeter Brune    Notes:
156162e0bf5SPeter Brune    The preferred calling sequence for users is to use SNESGetLineSearch() to acquire the SNESLineSearch instance
157162e0bf5SPeter Brune    already associated with the SNES.  This function is for developer use.
158f40b411bSPeter Brune 
159db781477SPatrick Sanan .seealso: `LineSearchDestroy()`, `SNESGetLineSearch()`
160f40b411bSPeter Brune @*/
161f40b411bSPeter Brune 
1629371c9d4SSatish Balay PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch) {
163f1c6b773SPeter Brune   SNESLineSearch linesearch;
164bf388a1fSBarry Smith 
165bf7f4e0aSPeter Brune   PetscFunctionBegin;
166ea5d4fccSPeter Brune   PetscValidPointer(outlinesearch, 2);
1679566063dSJacob Faibussowitsch   PetscCall(SNESInitializePackage());
1680298fd71SBarry Smith   *outlinesearch = NULL;
169f5af7f23SKarl Rupp 
1709566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(linesearch, SNESLINESEARCH_CLASSID, "SNESLineSearch", "Linesearch", "SNESLineSearch", comm, SNESLineSearchDestroy, SNESLineSearchView));
171bf7f4e0aSPeter Brune 
1720298fd71SBarry Smith   linesearch->vec_sol_new  = NULL;
1730298fd71SBarry Smith   linesearch->vec_func_new = NULL;
1740298fd71SBarry Smith   linesearch->vec_sol      = NULL;
1750298fd71SBarry Smith   linesearch->vec_func     = NULL;
1760298fd71SBarry Smith   linesearch->vec_update   = NULL;
1779bd66eb0SPeter Brune 
178bf7f4e0aSPeter Brune   linesearch->lambda       = 1.0;
179bf7f4e0aSPeter Brune   linesearch->fnorm        = 1.0;
180bf7f4e0aSPeter Brune   linesearch->ynorm        = 1.0;
181bf7f4e0aSPeter Brune   linesearch->xnorm        = 1.0;
182422a814eSBarry Smith   linesearch->result       = SNES_LINESEARCH_SUCCEEDED;
183bf7f4e0aSPeter Brune   linesearch->norms        = PETSC_TRUE;
184bf7f4e0aSPeter Brune   linesearch->keeplambda   = PETSC_FALSE;
185bf7f4e0aSPeter Brune   linesearch->damping      = 1.0;
186bf7f4e0aSPeter Brune   linesearch->maxstep      = 1e8;
187bf7f4e0aSPeter Brune   linesearch->steptol      = 1e-12;
188516fe3c3SPeter Brune   linesearch->rtol         = 1e-8;
189516fe3c3SPeter Brune   linesearch->atol         = 1e-15;
190516fe3c3SPeter Brune   linesearch->ltol         = 1e-8;
1910298fd71SBarry Smith   linesearch->precheckctx  = NULL;
1920298fd71SBarry Smith   linesearch->postcheckctx = NULL;
19359405d5eSPeter Brune   linesearch->max_its      = 1;
194bf7f4e0aSPeter Brune   linesearch->setupcalled  = PETSC_FALSE;
1953add74b1SBarry Smith   linesearch->monitor      = NULL;
196bf7f4e0aSPeter Brune   *outlinesearch           = linesearch;
197bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
198bf7f4e0aSPeter Brune }
199bf7f4e0aSPeter Brune 
200f40b411bSPeter Brune /*@
20178bcb3b5SPeter Brune    SNESLineSearchSetUp - Prepares the line search for being applied by allocating
20278bcb3b5SPeter Brune    any required vectors.
203f40b411bSPeter Brune 
204cd7522eaSPeter Brune    Collective on SNESLineSearch
205f40b411bSPeter Brune 
206f40b411bSPeter Brune    Input Parameters:
207f40b411bSPeter Brune .  linesearch - The LineSearch instance.
208f40b411bSPeter Brune 
209cd7522eaSPeter Brune    Notes:
210f190f2fcSBarry Smith    For most cases, this needn't be called by users or outside of SNESLineSearchApply().
211cd7522eaSPeter Brune    The only current case where this is called outside of this is for the VI
21278bcb3b5SPeter Brune    solvers, which modify the solution and work vectors before the first call
213cd7522eaSPeter Brune    of SNESLineSearchApply, requiring the SNESLineSearch work vectors to be
214cd7522eaSPeter Brune    allocated upfront.
215cd7522eaSPeter Brune 
21678bcb3b5SPeter Brune    Level: advanced
217f40b411bSPeter Brune 
218db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchReset()`
219f40b411bSPeter Brune @*/
220f40b411bSPeter Brune 
2219371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch) {
222bf7f4e0aSPeter Brune   PetscFunctionBegin;
223*48a46eb9SPierre Jolivet   if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC));
224bf7f4e0aSPeter Brune   if (!linesearch->setupcalled) {
225*48a46eb9SPierre Jolivet     if (!linesearch->vec_sol_new) PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new));
226*48a46eb9SPierre Jolivet     if (!linesearch->vec_func_new) PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new));
227dbbe0bcdSBarry Smith     if (linesearch->ops->setup) PetscUseTypeMethod(linesearch, setup);
2289566063dSJacob Faibussowitsch     if (!linesearch->ops->snesfunc) PetscCall(SNESLineSearchSetFunction(linesearch, SNESComputeFunction));
229bf7f4e0aSPeter Brune     linesearch->lambda      = linesearch->damping;
230bf7f4e0aSPeter Brune     linesearch->setupcalled = PETSC_TRUE;
231bf7f4e0aSPeter Brune   }
232bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
233bf7f4e0aSPeter Brune }
234bf7f4e0aSPeter Brune 
235f40b411bSPeter Brune /*@
236f190f2fcSBarry Smith    SNESLineSearchReset - Undoes the SNESLineSearchSetUp() and deletes any Vecs or Mats allocated by the line search.
237f40b411bSPeter Brune 
238f1c6b773SPeter Brune    Collective on SNESLineSearch
239f40b411bSPeter Brune 
240f40b411bSPeter Brune    Input Parameters:
241f40b411bSPeter Brune .  linesearch - The LineSearch instance.
242f40b411bSPeter Brune 
24395452b02SPatrick Sanan    Notes:
24495452b02SPatrick Sanan     Usually only called by SNESReset()
245f190f2fcSBarry Smith 
246f190f2fcSBarry Smith    Level: developer
247f40b411bSPeter Brune 
248db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchSetUp()`
249f40b411bSPeter Brune @*/
250f40b411bSPeter Brune 
2519371c9d4SSatish Balay PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch) {
252bf7f4e0aSPeter Brune   PetscFunctionBegin;
253dbbe0bcdSBarry Smith   if (linesearch->ops->reset) PetscUseTypeMethod(linesearch, reset);
254f5af7f23SKarl Rupp 
2559566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&linesearch->vec_sol_new));
2569566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&linesearch->vec_func_new));
257bf7f4e0aSPeter Brune 
2589566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(linesearch->nwork, &linesearch->work));
259f5af7f23SKarl Rupp 
260bf7f4e0aSPeter Brune   linesearch->nwork       = 0;
261bf7f4e0aSPeter Brune   linesearch->setupcalled = PETSC_FALSE;
262bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
263bf7f4e0aSPeter Brune }
264bf7f4e0aSPeter Brune 
265ed07d7d7SPeter Brune /*@C
266f190f2fcSBarry Smith    SNESLineSearchSetFunction - Sets the function evaluation used by the SNES line search
267ed07d7d7SPeter Brune 
268ed07d7d7SPeter Brune    Input Parameters:
269ed07d7d7SPeter Brune .  linesearch - the SNESLineSearch context
270f190f2fcSBarry Smith +  func       - function evaluation routine
271ed07d7d7SPeter Brune 
272ed07d7d7SPeter Brune    Level: developer
273ed07d7d7SPeter Brune 
27495452b02SPatrick Sanan    Notes:
27595452b02SPatrick Sanan     This is used internally by PETSc and not called by users
276f190f2fcSBarry Smith 
277db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESSetFunction()`
278ed07d7d7SPeter Brune @*/
2799371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES, Vec, Vec)) {
280ed07d7d7SPeter Brune   PetscFunctionBegin;
281ed07d7d7SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
282ed07d7d7SPeter Brune   linesearch->ops->snesfunc = func;
283ed07d7d7SPeter Brune   PetscFunctionReturn(0);
284ed07d7d7SPeter Brune }
285ed07d7d7SPeter Brune 
28686d74e61SPeter Brune /*@C
287f190f2fcSBarry Smith    SNESLineSearchSetPreCheck - Sets a user function that is called after the initial search direction has been computed but
288df3898eeSBarry Smith          before the line search routine has been applied. Allows the user to adjust the result of (usually a linear solve) that
289f190f2fcSBarry Smith          determined the search direction.
29086d74e61SPeter Brune 
291f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
29286d74e61SPeter Brune 
29386d74e61SPeter Brune    Input Parameters:
294f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
29584238204SBarry Smith .  func - [optional] function evaluation routine, see SNESLineSearchPreCheck() for the calling sequence
296f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
29786d74e61SPeter Brune 
29886d74e61SPeter Brune    Level: intermediate
29986d74e61SPeter Brune 
300f0b84518SBarry Smith    Notes:
301f0b84518SBarry Smith    Use `SNESLineSearchSetPostCheck()` to change the step after the line search.
302f0b84518SBarry Smith    search is complete.
303f0b84518SBarry Smith 
304f0b84518SBarry Smith    Use `SNESVISetVariableBounds()` and `SNESVISetComputeVariableBounds()` to cause `SNES` to automatically control the ranges of variables allowed.
305f0b84518SBarry Smith 
306f0b84518SBarry Smith .seealso: `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()`,
307f0b84518SBarry Smith           `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()
308f0b84518SBarry Smith 
30986d74e61SPeter Brune @*/
3109371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch, Vec, Vec, PetscBool *, void *), void *ctx) {
3119bd66eb0SPeter Brune   PetscFunctionBegin;
312f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
313f190f2fcSBarry Smith   if (func) linesearch->ops->precheck = func;
31486d74e61SPeter Brune   if (ctx) linesearch->precheckctx = ctx;
31586d74e61SPeter Brune   PetscFunctionReturn(0);
31686d74e61SPeter Brune }
31786d74e61SPeter Brune 
31886d74e61SPeter Brune /*@C
31953deb899SBarry Smith    SNESLineSearchGetPreCheck - Gets the pre-check function for the line search routine.
32086d74e61SPeter Brune 
321f899ff85SJose E. Roman    Input Parameter:
322f1c6b773SPeter Brune .  linesearch - the SNESLineSearch context
32386d74e61SPeter Brune 
32486d74e61SPeter Brune    Output Parameters:
32584238204SBarry Smith +  func       - [optional] function evaluation routine, see SNESLineSearchPreCheck() for calling sequence
326f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
32786d74e61SPeter Brune 
32886d74e61SPeter Brune    Level: intermediate
32986d74e61SPeter Brune 
330db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchSetPostCheck()`
33186d74e61SPeter Brune @*/
3329371c9d4SSatish Balay PetscErrorCode SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch, Vec, Vec, PetscBool *, void *), void **ctx) {
3339bd66eb0SPeter Brune   PetscFunctionBegin;
334f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
335f190f2fcSBarry Smith   if (func) *func = linesearch->ops->precheck;
33686d74e61SPeter Brune   if (ctx) *ctx = linesearch->precheckctx;
33786d74e61SPeter Brune   PetscFunctionReturn(0);
33886d74e61SPeter Brune }
33986d74e61SPeter Brune 
34086d74e61SPeter Brune /*@C
341f190f2fcSBarry Smith    SNESLineSearchSetPostCheck - Sets a user function that is called after the line search has been applied to determine the step
342f190f2fcSBarry Smith        direction and length. Allows the user a chance to change or override the decision of the line search routine
34386d74e61SPeter Brune 
344f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
34586d74e61SPeter Brune 
34686d74e61SPeter Brune    Input Parameters:
347f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
34884238204SBarry Smith .  func - [optional] function evaluation routine, see SNESLineSearchPostCheck()  for the calling sequence
349f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
35086d74e61SPeter Brune 
35186d74e61SPeter Brune    Level: intermediate
35286d74e61SPeter Brune 
353f0b84518SBarry Smith    Notes:
354f0b84518SBarry Smith    Use `SNESLineSearchSetPreCheck()` to change the step before the line search.
355f0b84518SBarry Smith    search is complete.
356f0b84518SBarry Smith 
357f0b84518SBarry Smith    Use `SNESVISetVariableBounds()` and `SNESVISetComputeVariableBounds()` to cause `SNES` to automatically control the ranges of variables allowed.
358f0b84518SBarry Smith 
359f0b84518SBarry Smith .seealso: `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchGetPostCheck()`,
360f0b84518SBarry Smith           `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()
36186d74e61SPeter Brune @*/
3629371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx) {
36386d74e61SPeter Brune   PetscFunctionBegin;
364f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
365f190f2fcSBarry Smith   if (func) linesearch->ops->postcheck = func;
36686d74e61SPeter Brune   if (ctx) linesearch->postcheckctx = ctx;
36786d74e61SPeter Brune   PetscFunctionReturn(0);
36886d74e61SPeter Brune }
36986d74e61SPeter Brune 
37086d74e61SPeter Brune /*@C
371f1c6b773SPeter Brune    SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine.
37286d74e61SPeter Brune 
373f899ff85SJose E. Roman    Input Parameter:
374f1c6b773SPeter Brune .  linesearch - the SNESLineSearch context
37586d74e61SPeter Brune 
37686d74e61SPeter Brune    Output Parameters:
37784238204SBarry Smith +  func - [optional] function evaluation routine, see for the calling sequence SNESLineSearchPostCheck()
378f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
37986d74e61SPeter Brune 
38086d74e61SPeter Brune    Level: intermediate
38186d74e61SPeter Brune 
382db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`
38386d74e61SPeter Brune @*/
3849371c9d4SSatish Balay PetscErrorCode SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx) {
3859bd66eb0SPeter Brune   PetscFunctionBegin;
386f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
387f190f2fcSBarry Smith   if (func) *func = linesearch->ops->postcheck;
38886d74e61SPeter Brune   if (ctx) *ctx = linesearch->postcheckctx;
38986d74e61SPeter Brune   PetscFunctionReturn(0);
39086d74e61SPeter Brune }
39186d74e61SPeter Brune 
392f40b411bSPeter Brune /*@
393f1c6b773SPeter Brune    SNESLineSearchPreCheck - Prepares the line search for being applied.
394f40b411bSPeter Brune 
395cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
396f40b411bSPeter Brune 
397f40b411bSPeter Brune    Input Parameters:
3987b1df9c1SPeter Brune +  linesearch - The linesearch instance.
3997b1df9c1SPeter Brune .  X - The current solution
4007b1df9c1SPeter Brune -  Y - The step direction
401f40b411bSPeter Brune 
402f40b411bSPeter Brune    Output Parameters:
4038e557f58SPeter Brune .  changed - Indicator that the precheck routine has changed anything
404f40b411bSPeter Brune 
405f0b84518SBarry Smith    Level: advanced
406f40b411bSPeter Brune 
407f0b84518SBarry Smith .seealso: `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`,
408f0b84518SBarry Smith           `SNESLineSearchGetPostCheck()``
409f40b411bSPeter Brune @*/
4109371c9d4SSatish Balay PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed) {
411bf7f4e0aSPeter Brune   PetscFunctionBegin;
412bf7f4e0aSPeter Brune   *changed = PETSC_FALSE;
4136b2b7091SBarry Smith   if (linesearch->ops->precheck) {
414dbbe0bcdSBarry Smith     PetscUseTypeMethod(linesearch, precheck, X, Y, changed, linesearch->precheckctx);
41538bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch, *changed, 4);
416bf7f4e0aSPeter Brune   }
417bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
418bf7f4e0aSPeter Brune }
419bf7f4e0aSPeter Brune 
420f40b411bSPeter Brune /*@
421ef46b1a6SStefano Zampini    SNESLineSearchPostCheck - Hook to modify step direction or updated solution after a successful linesearch
422f40b411bSPeter Brune 
423cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
424f40b411bSPeter Brune 
425f40b411bSPeter Brune    Input Parameters:
4267b1df9c1SPeter Brune +  linesearch - The linesearch context
4277b1df9c1SPeter Brune .  X - The last solution
4287b1df9c1SPeter Brune .  Y - The step direction
4297b1df9c1SPeter Brune -  W - The updated solution, W = X + lambda*Y for some lambda
430f40b411bSPeter Brune 
431f40b411bSPeter Brune    Output Parameters:
43278bcb3b5SPeter Brune +  changed_Y - Indicator if the direction Y has been changed.
43378bcb3b5SPeter Brune -  changed_W - Indicator if the new candidate solution W has been changed.
434f40b411bSPeter Brune 
435f190f2fcSBarry Smith    Level: developer
436f40b411bSPeter Brune 
437db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPrecheck()`, `SNESLineSearchGetPrecheck()`
438f40b411bSPeter Brune @*/
4399371c9d4SSatish Balay PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W) {
440bf7f4e0aSPeter Brune   PetscFunctionBegin;
441bf7f4e0aSPeter Brune   *changed_Y = PETSC_FALSE;
442bf7f4e0aSPeter Brune   *changed_W = PETSC_FALSE;
4436b2b7091SBarry Smith   if (linesearch->ops->postcheck) {
444dbbe0bcdSBarry Smith     PetscUseTypeMethod(linesearch, postcheck, X, Y, W, changed_Y, changed_W, linesearch->postcheckctx);
44538bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch, *changed_Y, 5);
44638bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch, *changed_W, 6);
44786d74e61SPeter Brune   }
44886d74e61SPeter Brune   PetscFunctionReturn(0);
44986d74e61SPeter Brune }
45086d74e61SPeter Brune 
45186d74e61SPeter Brune /*@C
45286d74e61SPeter Brune    SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration
45386d74e61SPeter Brune 
454cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
45586d74e61SPeter Brune 
4564165533cSJose E. Roman    Input Parameters:
45786d74e61SPeter Brune +  linesearch - linesearch context
45886d74e61SPeter Brune .  X - base state for this step
459907376e6SBarry Smith -  ctx - context for this function
46086d74e61SPeter Brune 
46197bb3fdcSJose E. Roman    Input/Output Parameter:
46297bb3fdcSJose E. Roman .  Y - correction, possibly modified
46397bb3fdcSJose E. Roman 
46497bb3fdcSJose E. Roman    Output Parameter:
46597bb3fdcSJose E. Roman .  changed - flag indicating that Y was modified
46686d74e61SPeter Brune 
46786d74e61SPeter Brune    Options Database Key:
468cd7522eaSPeter Brune +  -snes_linesearch_precheck_picard - activate this routine
469cd7522eaSPeter Brune -  -snes_linesearch_precheck_picard_angle - angle
47086d74e61SPeter Brune 
47186d74e61SPeter Brune    Level: advanced
47286d74e61SPeter Brune 
47386d74e61SPeter Brune    Notes:
47486d74e61SPeter Brune    This function should be passed to SNESLineSearchSetPreCheck()
47586d74e61SPeter Brune 
47686d74e61SPeter Brune    The justification for this method involves the linear convergence of a Picard iteration
47786d74e61SPeter Brune    so the Picard linearization should be provided in place of the "Jacobian". This correction
47886d74e61SPeter Brune    is generally not useful when using a Newton linearization.
47986d74e61SPeter Brune 
48086d74e61SPeter Brune    Reference:
48186d74e61SPeter Brune    Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology.
48286d74e61SPeter Brune 
483db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchSetPreCheck()`
48486d74e61SPeter Brune @*/
4859371c9d4SSatish Balay PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed, void *ctx) {
48686d74e61SPeter Brune   PetscReal   angle = *(PetscReal *)linesearch->precheckctx;
48786d74e61SPeter Brune   Vec         Ylast;
48886d74e61SPeter Brune   PetscScalar dot;
48986d74e61SPeter Brune   PetscInt    iter;
49086d74e61SPeter Brune   PetscReal   ynorm, ylastnorm, theta, angle_radians;
49186d74e61SPeter Brune   SNES        snes;
49286d74e61SPeter Brune 
49386d74e61SPeter Brune   PetscFunctionBegin;
4949566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
4959566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)snes, "SNESLineSearchPreCheckPicard_Ylast", (PetscObject *)&Ylast));
49686d74e61SPeter Brune   if (!Ylast) {
4979566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(Y, &Ylast));
4989566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)snes, "SNESLineSearchPreCheckPicard_Ylast", (PetscObject)Ylast));
4999566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)Ylast));
50086d74e61SPeter Brune   }
5019566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &iter));
50286d74e61SPeter Brune   if (iter < 2) {
5039566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y, Ylast));
50486d74e61SPeter Brune     *changed = PETSC_FALSE;
50586d74e61SPeter Brune     PetscFunctionReturn(0);
50686d74e61SPeter Brune   }
50786d74e61SPeter Brune 
5089566063dSJacob Faibussowitsch   PetscCall(VecDot(Y, Ylast, &dot));
5099566063dSJacob Faibussowitsch   PetscCall(VecNorm(Y, NORM_2, &ynorm));
5109566063dSJacob Faibussowitsch   PetscCall(VecNorm(Ylast, NORM_2, &ylastnorm));
51186d74e61SPeter Brune   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
512255453a1SBarry Smith   theta         = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm), -1.0, 1.0));
51386d74e61SPeter Brune   angle_radians = angle * PETSC_PI / 180.;
51486d74e61SPeter Brune   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
51586d74e61SPeter Brune     /* Modify the step Y */
51686d74e61SPeter Brune     PetscReal alpha, ydiffnorm;
5179566063dSJacob Faibussowitsch     PetscCall(VecAXPY(Ylast, -1.0, Y));
5189566063dSJacob Faibussowitsch     PetscCall(VecNorm(Ylast, NORM_2, &ydiffnorm));
519f85e2ce2SBarry Smith     alpha = (ydiffnorm > .001 * ylastnorm) ? ylastnorm / ydiffnorm : 1000.0;
5209566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y, Ylast));
5219566063dSJacob Faibussowitsch     PetscCall(VecScale(Y, alpha));
5229566063dSJacob 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));
523a47ec85fSBarry Smith     *changed = PETSC_TRUE;
52486d74e61SPeter Brune   } else {
5259566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Angle %14.12e degrees exceeds threshold %14.12e, no correction applied\n", (double)(theta * 180. / PETSC_PI), (double)angle));
5269566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y, Ylast));
52786d74e61SPeter Brune     *changed = PETSC_FALSE;
528bf7f4e0aSPeter Brune   }
529bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
530bf7f4e0aSPeter Brune }
531bf7f4e0aSPeter Brune 
532f40b411bSPeter Brune /*@
533cd7522eaSPeter Brune    SNESLineSearchApply - Computes the line-search update.
534f40b411bSPeter Brune 
535f1c6b773SPeter Brune    Collective on SNESLineSearch
536f40b411bSPeter Brune 
537f40b411bSPeter Brune    Input Parameters:
5388e557f58SPeter Brune +  linesearch - The linesearch context
5398e557f58SPeter Brune -  Y - The search direction
540f40b411bSPeter Brune 
5416b867d5aSJose E. Roman    Input/Output Parameters:
5426b867d5aSJose E. Roman +  X - The current solution, on output the new solution
5436b867d5aSJose E. Roman .  F - The current function, on output the new function
5446b867d5aSJose E. Roman -  fnorm - The current norm, on output the new function norm
545f40b411bSPeter Brune 
546cd7522eaSPeter Brune    Options Database Keys:
5470b00b554SBarry Smith + -snes_linesearch_type - basic (or equivalently none), bt, l2, cp, nleqerr, shell
548dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches
5491fe24845SBarry Smith . -snes_linesearch_damping - The linesearch damping parameter, default is 1.0 (no damping)
5501fe24845SBarry Smith . -snes_linesearch_norms   - Turn on/off the linesearch norms computation (SNESLineSearchSetComputeNorms())
5513c7d6663SPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess
5523c7d6663SPeter Brune - -snes_linesearch_max_it - The number of iterations for iterative line searches
553cd7522eaSPeter Brune 
554cd7522eaSPeter Brune    Notes:
555cd7522eaSPeter Brune    This is typically called from within a SNESSolve() implementation in order to
556cd7522eaSPeter Brune    help with convergence of the nonlinear method.  Various SNES types use line searches
557cd7522eaSPeter Brune    in different ways, but the overarching theme is that a line search is used to determine
558cd7522eaSPeter Brune    an optimal damping parameter of a step at each iteration of the method.  Each
55984238204SBarry Smith    application of the line search may invoke SNESComputeFunction() several times, and
560cd7522eaSPeter Brune    therefore may be fairly expensive.
561cd7522eaSPeter Brune 
562f40b411bSPeter Brune    Level: Intermediate
563f40b411bSPeter Brune 
564db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchPreCheck()`, `SNESLineSearchPostCheck()`, `SNESSolve()`, `SNESComputeFunction()`, `SNESLineSearchSetComputeNorms()`,
565db781477SPatrick Sanan           `SNESLineSearchType`, `SNESLineSearchSetType()`
566f40b411bSPeter Brune @*/
5679371c9d4SSatish Balay PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal *fnorm, Vec Y) {
568bf388a1fSBarry Smith   PetscFunctionBegin;
569f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
570bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
571bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
572064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(Y, VEC_CLASSID, 5);
573bf7f4e0aSPeter Brune 
574422a814eSBarry Smith   linesearch->result = SNES_LINESEARCH_SUCCEEDED;
575bf7f4e0aSPeter Brune 
576bf7f4e0aSPeter Brune   linesearch->vec_sol    = X;
577bf7f4e0aSPeter Brune   linesearch->vec_update = Y;
578bf7f4e0aSPeter Brune   linesearch->vec_func   = F;
579bf7f4e0aSPeter Brune 
5809566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchSetUp(linesearch));
581bf7f4e0aSPeter Brune 
582f5af7f23SKarl Rupp   if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */
583bf7f4e0aSPeter Brune 
584f5af7f23SKarl Rupp   if (fnorm) linesearch->fnorm = *fnorm;
585*48a46eb9SPierre Jolivet   else PetscCall(VecNorm(F, NORM_2, &linesearch->fnorm));
586bf7f4e0aSPeter Brune 
5879566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(SNESLINESEARCH_Apply, linesearch, X, F, Y));
588bf7f4e0aSPeter Brune 
589dbbe0bcdSBarry Smith   PetscUseTypeMethod(linesearch, apply);
590bf7f4e0aSPeter Brune 
5919566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(SNESLINESEARCH_Apply, linesearch, X, F, Y));
592bf7f4e0aSPeter Brune 
593f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
594bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
595bf7f4e0aSPeter Brune }
596bf7f4e0aSPeter Brune 
597f40b411bSPeter Brune /*@
598f1c6b773SPeter Brune    SNESLineSearchDestroy - Destroys the line search instance.
599f40b411bSPeter Brune 
600f1c6b773SPeter Brune    Collective on SNESLineSearch
601f40b411bSPeter Brune 
602f40b411bSPeter Brune    Input Parameters:
6038e557f58SPeter Brune .  linesearch - The linesearch context
604f40b411bSPeter Brune 
60584238204SBarry Smith    Level: developer
606f40b411bSPeter Brune 
607db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchReset()`, `SNESDestroy()`
608f40b411bSPeter Brune @*/
6099371c9d4SSatish Balay PetscErrorCode SNESLineSearchDestroy(SNESLineSearch *linesearch) {
610bf7f4e0aSPeter Brune   PetscFunctionBegin;
611bf7f4e0aSPeter Brune   if (!*linesearch) PetscFunctionReturn(0);
612f1c6b773SPeter Brune   PetscValidHeaderSpecific((*linesearch), SNESLINESEARCH_CLASSID, 1);
6139371c9d4SSatish Balay   if (--((PetscObject)(*linesearch))->refct > 0) {
6149371c9d4SSatish Balay     *linesearch = NULL;
6159371c9d4SSatish Balay     PetscFunctionReturn(0);
6169371c9d4SSatish Balay   }
6179566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsViewOff((PetscObject)*linesearch));
6189566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchReset(*linesearch));
619f5af7f23SKarl Rupp   if ((*linesearch)->ops->destroy) (*linesearch)->ops->destroy(*linesearch);
6209566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&(*linesearch)->monitor));
6219566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchMonitorCancel((*linesearch)));
6229566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(linesearch));
623bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
624bf7f4e0aSPeter Brune }
625bf7f4e0aSPeter Brune 
626f40b411bSPeter Brune /*@
627dcf2fd19SBarry Smith    SNESLineSearchSetDefaultMonitor - Turns on/off printing useful information and debugging output about the line search.
628bf7f4e0aSPeter Brune 
629bf7f4e0aSPeter Brune    Input Parameters:
630dcf2fd19SBarry Smith +  linesearch - the linesearch object
631dcf2fd19SBarry Smith -  viewer - an ASCII PetscViewer or NULL to turn off monitor
632bf7f4e0aSPeter Brune 
633dcf2fd19SBarry Smith    Logically Collective on SNESLineSearch
634bf7f4e0aSPeter Brune 
635bf7f4e0aSPeter Brune    Options Database:
636dcf2fd19SBarry Smith .   -snes_linesearch_monitor [:filename] - enables the monitor
637bf7f4e0aSPeter Brune 
638bf7f4e0aSPeter Brune    Level: intermediate
639bf7f4e0aSPeter Brune 
640dcf2fd19SBarry Smith    Developer Note: This monitor is implemented differently than the other SNESLineSearchMonitors that are set with
641d12e167eSBarry Smith      SNESLineSearchMonitorSet() since it is called in many locations of the line search routines to display aspects of the
642d12e167eSBarry Smith      line search that are not visible to the other monitors.
643dcf2fd19SBarry Smith 
644db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchGetDefaultMonitor()`, `PetscViewer`, `SNESLineSearchSetMonitor()`
645bf7f4e0aSPeter Brune @*/
6469371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetDefaultMonitor(SNESLineSearch linesearch, PetscViewer viewer) {
647bf7f4e0aSPeter Brune   PetscFunctionBegin;
6489566063dSJacob Faibussowitsch   if (viewer) PetscCall(PetscObjectReference((PetscObject)viewer));
6499566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&linesearch->monitor));
650dcf2fd19SBarry Smith   linesearch->monitor = viewer;
651bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
652bf7f4e0aSPeter Brune }
653bf7f4e0aSPeter Brune 
654f40b411bSPeter Brune /*@
655dcf2fd19SBarry Smith    SNESLineSearchGetDefaultMonitor - Gets the PetscViewer instance for the line search monitor.
6566a388c36SPeter Brune 
657f190f2fcSBarry Smith    Input Parameter:
6588e557f58SPeter Brune .  linesearch - linesearch context
659f40b411bSPeter Brune 
660f190f2fcSBarry Smith    Output Parameter:
6618e557f58SPeter Brune .  monitor - monitor context
662f40b411bSPeter Brune 
663f40b411bSPeter Brune    Logically Collective on SNES
664f40b411bSPeter Brune 
665f40b411bSPeter Brune    Options Database Keys:
6668e557f58SPeter Brune .   -snes_linesearch_monitor - enables the monitor
667f40b411bSPeter Brune 
668f40b411bSPeter Brune    Level: intermediate
669f40b411bSPeter Brune 
670db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchSetDefaultMonitor()`, `PetscViewer`
671f40b411bSPeter Brune @*/
6729371c9d4SSatish Balay PetscErrorCode SNESLineSearchGetDefaultMonitor(SNESLineSearch linesearch, PetscViewer *monitor) {
6736a388c36SPeter Brune   PetscFunctionBegin;
674f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
6756a388c36SPeter Brune   *monitor = linesearch->monitor;
6766a388c36SPeter Brune   PetscFunctionReturn(0);
6776a388c36SPeter Brune }
6786a388c36SPeter Brune 
679dcf2fd19SBarry Smith /*@C
680dcf2fd19SBarry Smith    SNESLineSearchMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
681dcf2fd19SBarry Smith 
682dcf2fd19SBarry Smith    Collective on SNESLineSearch
683dcf2fd19SBarry Smith 
684dcf2fd19SBarry Smith    Input Parameters:
685dcf2fd19SBarry Smith +  ls - LineSearch object you wish to monitor
686dcf2fd19SBarry Smith .  name - the monitor type one is seeking
687dcf2fd19SBarry Smith .  help - message indicating what monitoring is done
688dcf2fd19SBarry Smith .  manual - manual page for the monitor
689dcf2fd19SBarry Smith .  monitor - the monitor function
690dcf2fd19SBarry 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 objects
691dcf2fd19SBarry Smith 
692dcf2fd19SBarry Smith    Level: developer
693dcf2fd19SBarry Smith 
694db781477SPatrick Sanan .seealso: `PetscOptionsGetViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
695db781477SPatrick Sanan           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
696db781477SPatrick Sanan           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
697db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
698c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
699db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
700db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
701dcf2fd19SBarry Smith @*/
7029371c9d4SSatish Balay PetscErrorCode SNESLineSearchMonitorSetFromOptions(SNESLineSearch ls, const char name[], const char help[], const char manual[], PetscErrorCode (*monitor)(SNESLineSearch, PetscViewerAndFormat *), PetscErrorCode (*monitorsetup)(SNESLineSearch, PetscViewerAndFormat *)) {
703dcf2fd19SBarry Smith   PetscViewer       viewer;
704dcf2fd19SBarry Smith   PetscViewerFormat format;
705dcf2fd19SBarry Smith   PetscBool         flg;
706dcf2fd19SBarry Smith 
707dcf2fd19SBarry Smith   PetscFunctionBegin;
7089566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)ls), ((PetscObject)ls)->options, ((PetscObject)ls)->prefix, name, &viewer, &format, &flg));
709dcf2fd19SBarry Smith   if (flg) {
710d12e167eSBarry Smith     PetscViewerAndFormat *vf;
7119566063dSJacob Faibussowitsch     PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf));
7129566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)viewer));
7131baa6e33SBarry Smith     if (monitorsetup) PetscCall((*monitorsetup)(ls, vf));
7149566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchMonitorSet(ls, (PetscErrorCode(*)(SNESLineSearch, void *))monitor, vf, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy));
715dcf2fd19SBarry Smith   }
716dcf2fd19SBarry Smith   PetscFunctionReturn(0);
717dcf2fd19SBarry Smith }
718dcf2fd19SBarry Smith 
719f40b411bSPeter Brune /*@
720f1c6b773SPeter Brune    SNESLineSearchSetFromOptions - Sets options for the line search
721f40b411bSPeter Brune 
722f899ff85SJose E. Roman    Input Parameter:
7238e557f58SPeter Brune .  linesearch - linesearch context
724f40b411bSPeter Brune 
725f40b411bSPeter Brune    Options Database Keys:
7260b00b554SBarry Smith + -snes_linesearch_type <type> - basic (or equivalently none), bt, l2, cp, nleqerr, shell
7273c7d6663SPeter Brune . -snes_linesearch_order <order> - 1, 2, 3.  Most types only support certain orders (bt supports 2 or 3)
7281fe24845SBarry Smith . -snes_linesearch_norms   - Turn on/off the linesearch norms for the basic linesearch typem (SNESLineSearchSetComputeNorms())
72971eef1aeSPeter Brune . -snes_linesearch_minlambda - The minimum step length
7301a9b3a06SPeter Brune . -snes_linesearch_maxstep - The maximum step size
7311a9b3a06SPeter Brune . -snes_linesearch_rtol - Relative tolerance for iterative line searches
7321a9b3a06SPeter Brune . -snes_linesearch_atol - Absolute tolerance for iterative line searches
7331a9b3a06SPeter Brune . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches
7341a9b3a06SPeter Brune . -snes_linesearch_max_it - The number of iterations for iterative line searches
735dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches
736dcf2fd19SBarry Smith . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine
7378e557f58SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter
738cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess.
7391a9b3a06SPeter Brune . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method
740d8d34be6SBarry Smith - -snes_linesearch_precheck_picard_angle - Angle used in Picard precheck method
741f40b411bSPeter Brune 
742f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
743f40b411bSPeter Brune 
744f40b411bSPeter Brune    Level: intermediate
745f40b411bSPeter Brune 
746db781477SPatrick Sanan .seealso: `SNESLineSearchCreate()`, `SNESLineSearchSetOrder()`, `SNESLineSearchSetType()`, `SNESLineSearchSetTolerances()`, `SNESLineSearchSetDamping()`, `SNESLineSearchPreCheckPicard()`,
747db781477SPatrick Sanan           `SNESLineSearchType`, `SNESLineSearchSetComputeNorms()`
748f40b411bSPeter Brune @*/
7499371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch) {
7501a4f838cSPeter Brune   const char *deft = SNESLINESEARCHBASIC;
751bf7f4e0aSPeter Brune   char        type[256];
752bf7f4e0aSPeter Brune   PetscBool   flg, set;
753dcf2fd19SBarry Smith   PetscViewer viewer;
754bf388a1fSBarry Smith 
755bf7f4e0aSPeter Brune   PetscFunctionBegin;
7569566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchRegisterAll());
757bf7f4e0aSPeter Brune 
758d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)linesearch);
759f5af7f23SKarl Rupp   if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name;
7609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-snes_linesearch_type", "Linesearch type", "SNESLineSearchSetType", SNESLineSearchList, deft, type, 256, &flg));
761bf7f4e0aSPeter Brune   if (flg) {
7629566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetType(linesearch, type));
763bf7f4e0aSPeter Brune   } else if (!((PetscObject)linesearch)->type_name) {
7649566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetType(linesearch, deft));
765bf7f4e0aSPeter Brune   }
766bf7f4e0aSPeter Brune 
7679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)linesearch), ((PetscObject)linesearch)->options, ((PetscObject)linesearch)->prefix, "-snes_linesearch_monitor", &viewer, NULL, &set));
768dcf2fd19SBarry Smith   if (set) {
7699566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetDefaultMonitor(linesearch, viewer));
7709566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
771dcf2fd19SBarry Smith   }
7729566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchMonitorSetFromOptions(linesearch, "-snes_linesearch_monitor_solution_update", "View correction at each iteration", "SNESLineSearchMonitorSolutionUpdate", SNESLineSearchMonitorSolutionUpdate, NULL));
773bf7f4e0aSPeter Brune 
7741a9b3a06SPeter Brune   /* tolerances */
7759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_minlambda", "Minimum step length", "SNESLineSearchSetTolerances", linesearch->steptol, &linesearch->steptol, NULL));
7769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_maxstep", "Maximum step size", "SNESLineSearchSetTolerances", linesearch->maxstep, &linesearch->maxstep, NULL));
7779566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_rtol", "Relative tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->rtol, &linesearch->rtol, NULL));
7789566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_atol", "Absolute tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->atol, &linesearch->atol, NULL));
7799566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_ltol", "Change in lambda tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->ltol, &linesearch->ltol, NULL));
7809566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_linesearch_max_it", "Maximum iterations for iterative line searches", "SNESLineSearchSetTolerances", linesearch->max_its, &linesearch->max_its, NULL));
7817a35526eSPeter Brune 
7821a9b3a06SPeter Brune   /* damping parameters */
7839566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_damping", "Line search damping and initial step guess", "SNESLineSearchSetDamping", linesearch->damping, &linesearch->damping, NULL));
7841a9b3a06SPeter Brune 
7859566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_linesearch_keeplambda", "Use previous lambda as damping", "SNESLineSearchSetKeepLambda", linesearch->keeplambda, &linesearch->keeplambda, NULL));
7861a9b3a06SPeter Brune 
7871a9b3a06SPeter Brune   /* precheck */
7889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_linesearch_precheck_picard", "Use a correction that sometimes improves convergence of Picard iteration", "SNESLineSearchPreCheckPicard", flg, &flg, &set));
7897a35526eSPeter Brune   if (set) {
7907a35526eSPeter Brune     if (flg) {
7917a35526eSPeter Brune       linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */
792f5af7f23SKarl Rupp 
793d0609cedSBarry 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));
7949566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchSetPreCheck(linesearch, SNESLineSearchPreCheckPicard, &linesearch->precheck_picard_angle));
7957a35526eSPeter Brune     } else {
7969566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchSetPreCheck(linesearch, NULL, NULL));
7977a35526eSPeter Brune     }
7987a35526eSPeter Brune   }
7999566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_linesearch_order", "Order of approximation used in the line search", "SNESLineSearchSetOrder", linesearch->order, &linesearch->order, NULL));
8009566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_linesearch_norms", "Compute final norms in line search", "SNESLineSearchSetComputeNorms", linesearch->norms, &linesearch->norms, NULL));
8017a35526eSPeter Brune 
802dbbe0bcdSBarry Smith   PetscTryTypeMethod(linesearch, setfromoptions, PetscOptionsObject);
8035a9b6599SPeter Brune 
804dbbe0bcdSBarry Smith   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)linesearch, PetscOptionsObject));
805d0609cedSBarry Smith   PetscOptionsEnd();
806bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
807bf7f4e0aSPeter Brune }
808bf7f4e0aSPeter Brune 
809f40b411bSPeter Brune /*@
810f190f2fcSBarry Smith    SNESLineSearchView - Prints useful information about the line search
811f40b411bSPeter Brune 
812f40b411bSPeter Brune    Input Parameters:
8138e557f58SPeter Brune .  linesearch - linesearch context
814f40b411bSPeter Brune 
815f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
816f40b411bSPeter Brune 
817f40b411bSPeter Brune    Level: intermediate
818f40b411bSPeter Brune 
819db781477SPatrick Sanan .seealso: `SNESLineSearchCreate()`
820f40b411bSPeter Brune @*/
8219371c9d4SSatish Balay PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer) {
8227f1410a3SPeter Brune   PetscBool iascii;
823bf388a1fSBarry Smith 
824bf7f4e0aSPeter Brune   PetscFunctionBegin;
8257f1410a3SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
826*48a46eb9SPierre Jolivet   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch), &viewer));
8277f1410a3SPeter Brune   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
8287f1410a3SPeter Brune   PetscCheckSameComm(linesearch, 1, viewer, 2);
829f40b411bSPeter Brune 
8309566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
8317f1410a3SPeter Brune   if (iascii) {
8329566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)linesearch, viewer));
8339566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
834dbbe0bcdSBarry Smith     PetscTryTypeMethod(linesearch, view, viewer);
8359566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
8369566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep, (double)linesearch->steptol));
8379566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol, (double)linesearch->atol, (double)linesearch->ltol));
83863a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  maximum iterations=%" PetscInt_FMT "\n", linesearch->max_its));
8396b2b7091SBarry Smith     if (linesearch->ops->precheck) {
8406b2b7091SBarry Smith       if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) {
84163a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  using precheck step to speed up Picard convergence\n"));
8427f1410a3SPeter Brune       } else {
84363a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  using user-defined precheck step\n"));
8447f1410a3SPeter Brune       }
8457f1410a3SPeter Brune     }
846*48a46eb9SPierre Jolivet     if (linesearch->ops->postcheck) PetscCall(PetscViewerASCIIPrintf(viewer, "  using user-defined postcheck step\n"));
8477f1410a3SPeter Brune   }
848bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
849bf7f4e0aSPeter Brune }
850bf7f4e0aSPeter Brune 
851ea5d4fccSPeter Brune /*@C
852a80ff896SJed Brown    SNESLineSearchGetType - Gets the linesearch type
853a80ff896SJed Brown 
854a80ff896SJed Brown    Logically Collective on SNESLineSearch
855a80ff896SJed Brown 
856a80ff896SJed Brown    Input Parameters:
857a80ff896SJed Brown .  linesearch - linesearch context
858a80ff896SJed Brown 
859a80ff896SJed Brown    Output Parameters:
860a80ff896SJed Brown -  type - The type of line search, or NULL if not set
861a80ff896SJed Brown 
862a80ff896SJed Brown    Level: intermediate
863a80ff896SJed Brown 
864db781477SPatrick Sanan .seealso: `SNESLineSearchCreate()`, `SNESLineSearchType`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchSetType()`
865a80ff896SJed Brown @*/
8669371c9d4SSatish Balay PetscErrorCode SNESLineSearchGetType(SNESLineSearch linesearch, SNESLineSearchType *type) {
867a80ff896SJed Brown   PetscFunctionBegin;
868a80ff896SJed Brown   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
869dadcf809SJacob Faibussowitsch   PetscValidPointer(type, 2);
870a80ff896SJed Brown   *type = ((PetscObject)linesearch)->type_name;
871a80ff896SJed Brown   PetscFunctionReturn(0);
872a80ff896SJed Brown }
873a80ff896SJed Brown 
874a80ff896SJed Brown /*@C
875f1c6b773SPeter Brune    SNESLineSearchSetType - Sets the linesearch type
876f40b411bSPeter Brune 
877f190f2fcSBarry Smith    Logically Collective on SNESLineSearch
878f190f2fcSBarry Smith 
879f40b411bSPeter Brune    Input Parameters:
8808e557f58SPeter Brune +  linesearch - linesearch context
881f40b411bSPeter Brune -  type - The type of line search to be used
882f40b411bSPeter Brune 
883cd7522eaSPeter Brune    Available Types:
8840b00b554SBarry Smith +  SNESLINESEARCHBASIC - (or equivalently SNESLINESEARCHNONE) Simple damping line search, defaults to using the full Newton step
8851fe24845SBarry Smith .  SNESLINESEARCHBT - Backtracking line search over the L2 norm of the function
8861fe24845SBarry Smith .  SNESLINESEARCHL2 - Secant line search over the L2 norm of the function
8871fe24845SBarry Smith .  SNESLINESEARCHCP - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x)
8881fe24845SBarry Smith .  SNESLINESEARCHNLEQERR - Affine-covariant error-oriented linesearch
8891fe24845SBarry Smith -  SNESLINESEARCHSHELL - User provided SNESLineSearch implementation
8901fe24845SBarry Smith 
8911fe24845SBarry Smith    Options Database:
8920b00b554SBarry Smith .  -snes_linesearch_type <type> - basic (or equivalently none), bt, l2, cp, nleqerr, shell
893cd7522eaSPeter Brune 
894f40b411bSPeter Brune    Level: intermediate
895f40b411bSPeter Brune 
896db781477SPatrick Sanan .seealso: `SNESLineSearchCreate()`, `SNESLineSearchType`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchGetType()`
897f40b411bSPeter Brune @*/
8989371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type) {
899bf7f4e0aSPeter Brune   PetscBool match;
9005f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(SNESLineSearch);
901bf7f4e0aSPeter Brune 
902bf7f4e0aSPeter Brune   PetscFunctionBegin;
903f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
904bf7f4e0aSPeter Brune   PetscValidCharPointer(type, 2);
905bf7f4e0aSPeter Brune 
9069566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)linesearch, type, &match));
907bf7f4e0aSPeter Brune   if (match) PetscFunctionReturn(0);
908bf7f4e0aSPeter Brune 
9099566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(SNESLineSearchList, type, &r));
9105f80ce2aSJacob Faibussowitsch   PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested Line Search type %s", type);
911bf7f4e0aSPeter Brune   /* Destroy the previous private linesearch context */
912bf7f4e0aSPeter Brune   if (linesearch->ops->destroy) {
9139566063dSJacob Faibussowitsch     PetscCall((*(linesearch)->ops->destroy)(linesearch));
9140298fd71SBarry Smith     linesearch->ops->destroy = NULL;
915bf7f4e0aSPeter Brune   }
916f1c6b773SPeter Brune   /* Reinitialize function pointers in SNESLineSearchOps structure */
9179e5d0892SLisandro Dalcin   linesearch->ops->apply          = NULL;
9189e5d0892SLisandro Dalcin   linesearch->ops->view           = NULL;
9199e5d0892SLisandro Dalcin   linesearch->ops->setfromoptions = NULL;
9209e5d0892SLisandro Dalcin   linesearch->ops->destroy        = NULL;
921bf7f4e0aSPeter Brune 
9229566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)linesearch, type));
9239566063dSJacob Faibussowitsch   PetscCall((*r)(linesearch));
924bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
925bf7f4e0aSPeter Brune }
926bf7f4e0aSPeter Brune 
927f40b411bSPeter Brune /*@
92878bcb3b5SPeter Brune    SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation.
929f40b411bSPeter Brune 
930f40b411bSPeter Brune    Input Parameters:
9318e557f58SPeter Brune +  linesearch - linesearch context
932f40b411bSPeter Brune -  snes - The snes instance
933f40b411bSPeter Brune 
93478bcb3b5SPeter Brune    Level: developer
93578bcb3b5SPeter Brune 
93678bcb3b5SPeter Brune    Notes:
937f190f2fcSBarry Smith    This happens automatically when the line search is obtained/created with
9387601faf0SJed Brown    SNESGetLineSearch().  This routine is therefore mainly called within SNES
93978bcb3b5SPeter Brune    implementations.
940f40b411bSPeter Brune 
9418141a3b9SPeter Brune    Level: developer
942f40b411bSPeter Brune 
943db781477SPatrick Sanan .seealso: `SNESLineSearchGetSNES()`, `SNESLineSearchSetVecs()`, `SNES`
944f40b411bSPeter Brune @*/
9459371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes) {
946bf7f4e0aSPeter Brune   PetscFunctionBegin;
947f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
948bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 2);
949bf7f4e0aSPeter Brune   linesearch->snes = snes;
950bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
951bf7f4e0aSPeter Brune }
952bf7f4e0aSPeter Brune 
953f40b411bSPeter Brune /*@
9548141a3b9SPeter Brune    SNESLineSearchGetSNES - Gets the SNES instance associated with the line search.
9558141a3b9SPeter Brune    Having an associated SNES is necessary because most line search implementations must be able to
9568141a3b9SPeter Brune    evaluate the function using SNESComputeFunction() for the associated SNES.  This routine
9578141a3b9SPeter Brune    is used in the line search implementations when one must get this associated SNES instance.
958f40b411bSPeter Brune 
959f40b411bSPeter Brune    Input Parameters:
9608e557f58SPeter Brune .  linesearch - linesearch context
961f40b411bSPeter Brune 
962f40b411bSPeter Brune    Output Parameters:
963f40b411bSPeter Brune .  snes - The snes instance
964f40b411bSPeter Brune 
9658141a3b9SPeter Brune    Level: developer
966f40b411bSPeter Brune 
967db781477SPatrick Sanan .seealso: `SNESLineSearchGetSNES()`, `SNESLineSearchSetVecs()`, `SNES`
968f40b411bSPeter Brune @*/
9699371c9d4SSatish Balay PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes) {
970bf7f4e0aSPeter Brune   PetscFunctionBegin;
971f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
9726a388c36SPeter Brune   PetscValidPointer(snes, 2);
973bf7f4e0aSPeter Brune   *snes = linesearch->snes;
974bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
975bf7f4e0aSPeter Brune }
976bf7f4e0aSPeter Brune 
977f40b411bSPeter Brune /*@
978f1c6b773SPeter Brune    SNESLineSearchGetLambda - Gets the last linesearch steplength discovered.
979f40b411bSPeter Brune 
980f40b411bSPeter Brune    Input Parameters:
9818e557f58SPeter Brune .  linesearch - linesearch context
982f40b411bSPeter Brune 
983f40b411bSPeter Brune    Output Parameters:
984cd7522eaSPeter Brune .  lambda - The last steplength computed during SNESLineSearchApply()
985f40b411bSPeter Brune 
98678bcb3b5SPeter Brune    Level: advanced
98778bcb3b5SPeter Brune 
9888e557f58SPeter Brune    Notes:
9898e557f58SPeter Brune    This is useful in methods where the solver is ill-scaled and
99078bcb3b5SPeter Brune    requires some adaptive notion of the difference in scale between the
99178bcb3b5SPeter Brune    solution and the function.  For instance, SNESQN may be scaled by the
99278bcb3b5SPeter Brune    line search lambda using the argument -snes_qn_scaling ls.
99378bcb3b5SPeter Brune 
994db781477SPatrick Sanan .seealso: `SNESLineSearchSetLambda()`, `SNESLineSearchGetDamping()`, `SNESLineSearchApply()`
995f40b411bSPeter Brune @*/
9969371c9d4SSatish Balay PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch linesearch, PetscReal *lambda) {
9976a388c36SPeter Brune   PetscFunctionBegin;
998f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
999534a8f05SLisandro Dalcin   PetscValidRealPointer(lambda, 2);
10006a388c36SPeter Brune   *lambda = linesearch->lambda;
10016a388c36SPeter Brune   PetscFunctionReturn(0);
10026a388c36SPeter Brune }
10036a388c36SPeter Brune 
1004f40b411bSPeter Brune /*@
1005f1c6b773SPeter Brune    SNESLineSearchSetLambda - Sets the linesearch steplength.
1006f40b411bSPeter Brune 
1007f40b411bSPeter Brune    Input Parameters:
10088e557f58SPeter Brune +  linesearch - linesearch context
1009f40b411bSPeter Brune -  lambda - The last steplength.
1010f40b411bSPeter Brune 
1011cd7522eaSPeter Brune    Notes:
1012f190f2fcSBarry Smith    This routine is typically used within implementations of SNESLineSearchApply()
1013cd7522eaSPeter Brune    to set the final steplength.  This routine (and SNESLineSearchGetLambda()) were
1014cd7522eaSPeter Brune    added in order to facilitate Quasi-Newton methods that use the previous steplength
1015cd7522eaSPeter Brune    as an inner scaling parameter.
1016cd7522eaSPeter Brune 
101778bcb3b5SPeter Brune    Level: advanced
1018f40b411bSPeter Brune 
1019db781477SPatrick Sanan .seealso: `SNESLineSearchGetLambda()`
1020f40b411bSPeter Brune @*/
10219371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda) {
10226a388c36SPeter Brune   PetscFunctionBegin;
1023f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
10246a388c36SPeter Brune   linesearch->lambda = lambda;
10256a388c36SPeter Brune   PetscFunctionReturn(0);
10266a388c36SPeter Brune }
10276a388c36SPeter Brune 
1028f40b411bSPeter Brune /*@
10293c7d6663SPeter Brune    SNESLineSearchGetTolerances - Gets the tolerances for the linesearch.  These include
103078bcb3b5SPeter Brune    tolerances for the relative and absolute change in the function norm, the change
103178bcb3b5SPeter Brune    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
103278bcb3b5SPeter Brune    and the maximum number of iterations the line search procedure may take.
1033f40b411bSPeter Brune 
1034f899ff85SJose E. Roman    Input Parameter:
10358e557f58SPeter Brune .  linesearch - linesearch context
1036f40b411bSPeter Brune 
1037f40b411bSPeter Brune    Output Parameters:
1038516fe3c3SPeter Brune +  steptol - The minimum steplength
10396cc8e53bSPeter Brune .  maxstep - The maximum steplength
1040516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
1041516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
1042516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
1043516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
1044f40b411bSPeter Brune 
104578bcb3b5SPeter Brune    Level: intermediate
104678bcb3b5SPeter Brune 
104778bcb3b5SPeter Brune    Notes:
104878bcb3b5SPeter Brune    Different line searches may implement these parameters slightly differently as
10493c7d6663SPeter Brune    the type requires.
1050516fe3c3SPeter Brune 
1051db781477SPatrick Sanan .seealso: `SNESLineSearchSetTolerances()`
1052f40b411bSPeter Brune @*/
10539371c9d4SSatish Balay PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch linesearch, PetscReal *steptol, PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its) {
10546a388c36SPeter Brune   PetscFunctionBegin;
1055f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1056516fe3c3SPeter Brune   if (steptol) {
1057534a8f05SLisandro Dalcin     PetscValidRealPointer(steptol, 2);
10586a388c36SPeter Brune     *steptol = linesearch->steptol;
1059516fe3c3SPeter Brune   }
1060516fe3c3SPeter Brune   if (maxstep) {
1061534a8f05SLisandro Dalcin     PetscValidRealPointer(maxstep, 3);
1062516fe3c3SPeter Brune     *maxstep = linesearch->maxstep;
1063516fe3c3SPeter Brune   }
1064516fe3c3SPeter Brune   if (rtol) {
1065534a8f05SLisandro Dalcin     PetscValidRealPointer(rtol, 4);
1066516fe3c3SPeter Brune     *rtol = linesearch->rtol;
1067516fe3c3SPeter Brune   }
1068516fe3c3SPeter Brune   if (atol) {
1069534a8f05SLisandro Dalcin     PetscValidRealPointer(atol, 5);
1070516fe3c3SPeter Brune     *atol = linesearch->atol;
1071516fe3c3SPeter Brune   }
1072516fe3c3SPeter Brune   if (ltol) {
1073534a8f05SLisandro Dalcin     PetscValidRealPointer(ltol, 6);
1074516fe3c3SPeter Brune     *ltol = linesearch->ltol;
1075516fe3c3SPeter Brune   }
1076516fe3c3SPeter Brune   if (max_its) {
1077534a8f05SLisandro Dalcin     PetscValidIntPointer(max_its, 7);
1078516fe3c3SPeter Brune     *max_its = linesearch->max_its;
1079516fe3c3SPeter Brune   }
10806a388c36SPeter Brune   PetscFunctionReturn(0);
10816a388c36SPeter Brune }
10826a388c36SPeter Brune 
1083f40b411bSPeter Brune /*@
10843c7d6663SPeter Brune    SNESLineSearchSetTolerances -  Gets the tolerances for the linesearch.  These include
108578bcb3b5SPeter Brune    tolerances for the relative and absolute change in the function norm, the change
108678bcb3b5SPeter Brune    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
108778bcb3b5SPeter Brune    and the maximum number of iterations the line search procedure may take.
1088f40b411bSPeter Brune 
1089f40b411bSPeter Brune    Input Parameters:
10908e557f58SPeter Brune +  linesearch - linesearch context
1091516fe3c3SPeter Brune .  steptol - The minimum steplength
10926cc8e53bSPeter Brune .  maxstep - The maximum steplength
1093516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
1094516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
1095516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
1096516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
1097f40b411bSPeter Brune 
109878bcb3b5SPeter Brune    Notes:
10993c7d6663SPeter Brune    The user may choose to not set any of the tolerances using PETSC_DEFAULT in
110078bcb3b5SPeter Brune    place of an argument.
1101f40b411bSPeter Brune 
110278bcb3b5SPeter Brune    Level: intermediate
1103516fe3c3SPeter Brune 
1104db781477SPatrick Sanan .seealso: `SNESLineSearchGetTolerances()`
1105f40b411bSPeter Brune @*/
11069371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch linesearch, PetscReal steptol, PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its) {
11076a388c36SPeter Brune   PetscFunctionBegin;
1108f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1109d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, steptol, 2);
1110d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, maxstep, 3);
1111d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, rtol, 4);
1112d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, atol, 5);
1113d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, ltol, 6);
1114d3952184SSatish Balay   PetscValidLogicalCollectiveInt(linesearch, max_its, 7);
1115d3952184SSatish Balay 
1116d3952184SSatish Balay   if (steptol != PETSC_DEFAULT) {
11175f80ce2aSJacob Faibussowitsch     PetscCheck(steptol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Minimum step length %14.12e must be non-negative", (double)steptol);
11186a388c36SPeter Brune     linesearch->steptol = steptol;
1119d3952184SSatish Balay   }
1120d3952184SSatish Balay 
1121d3952184SSatish Balay   if (maxstep != PETSC_DEFAULT) {
11225f80ce2aSJacob Faibussowitsch     PetscCheck(maxstep >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum step length %14.12e must be non-negative", (double)maxstep);
1123516fe3c3SPeter Brune     linesearch->maxstep = maxstep;
1124d3952184SSatish Balay   }
1125d3952184SSatish Balay 
1126d3952184SSatish Balay   if (rtol != PETSC_DEFAULT) {
11272061ca32SJunchao 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);
1128516fe3c3SPeter Brune     linesearch->rtol = rtol;
1129d3952184SSatish Balay   }
1130d3952184SSatish Balay 
1131d3952184SSatish Balay   if (atol != PETSC_DEFAULT) {
11325f80ce2aSJacob Faibussowitsch     PetscCheck(atol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Absolute tolerance %14.12e must be non-negative", (double)atol);
1133516fe3c3SPeter Brune     linesearch->atol = atol;
1134d3952184SSatish Balay   }
1135d3952184SSatish Balay 
1136d3952184SSatish Balay   if (ltol != PETSC_DEFAULT) {
11375f80ce2aSJacob Faibussowitsch     PetscCheck(ltol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Lambda tolerance %14.12e must be non-negative", (double)ltol);
1138516fe3c3SPeter Brune     linesearch->ltol = ltol;
1139d3952184SSatish Balay   }
1140d3952184SSatish Balay 
1141d3952184SSatish Balay   if (max_its != PETSC_DEFAULT) {
114263a3b9bcSJacob Faibussowitsch     PetscCheck(max_its >= 0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of iterations %" PetscInt_FMT " must be non-negative", max_its);
1143516fe3c3SPeter Brune     linesearch->max_its = max_its;
1144d3952184SSatish Balay   }
11456a388c36SPeter Brune   PetscFunctionReturn(0);
11466a388c36SPeter Brune }
11476a388c36SPeter Brune 
1148f40b411bSPeter Brune /*@
1149f1c6b773SPeter Brune    SNESLineSearchGetDamping - Gets the line search damping parameter.
1150f40b411bSPeter Brune 
1151f40b411bSPeter Brune    Input Parameters:
11528e557f58SPeter Brune .  linesearch - linesearch context
1153f40b411bSPeter Brune 
1154f40b411bSPeter Brune    Output Parameters:
11558e557f58SPeter Brune .  damping - The damping parameter
1156f40b411bSPeter Brune 
115778bcb3b5SPeter Brune    Level: advanced
1158f40b411bSPeter Brune 
1159db781477SPatrick Sanan .seealso: `SNESLineSearchGetStepTolerance()`, `SNESQN`
1160f40b411bSPeter Brune @*/
1161f40b411bSPeter Brune 
11629371c9d4SSatish Balay PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch linesearch, PetscReal *damping) {
11636a388c36SPeter Brune   PetscFunctionBegin;
1164f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1165534a8f05SLisandro Dalcin   PetscValidRealPointer(damping, 2);
11666a388c36SPeter Brune   *damping = linesearch->damping;
11676a388c36SPeter Brune   PetscFunctionReturn(0);
11686a388c36SPeter Brune }
11696a388c36SPeter Brune 
1170f40b411bSPeter Brune /*@
1171fd292e60Sprj-    SNESLineSearchSetDamping - Sets the line search damping parameter.
1172f40b411bSPeter Brune 
1173f40b411bSPeter Brune    Input Parameters:
117403fd4120SBarry Smith +  linesearch - linesearch context
117503fd4120SBarry Smith -  damping - The damping parameter
1176f40b411bSPeter Brune 
117703fd4120SBarry Smith    Options Database:
117803fd4120SBarry Smith .   -snes_linesearch_damping
1179f40b411bSPeter Brune    Level: intermediate
1180f40b411bSPeter Brune 
1181cd7522eaSPeter Brune    Notes:
11820b00b554SBarry Smith    The basic (also known as the none) line search merely takes the update step scaled by the damping parameter.
1183cd7522eaSPeter Brune    The use of the damping parameter in the l2 and cp line searches is much more subtle;
118478bcb3b5SPeter Brune    it is used as a starting point in calculating the secant step. However, the eventual
1185cd7522eaSPeter Brune    step may be of greater length than the damping parameter.  In the bt line search it is
1186cd7522eaSPeter Brune    used as the maximum possible step length, as the bt line search only backtracks.
1187cd7522eaSPeter Brune 
1188db781477SPatrick Sanan .seealso: `SNESLineSearchGetDamping()`
1189f40b411bSPeter Brune @*/
11909371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch linesearch, PetscReal damping) {
11916a388c36SPeter Brune   PetscFunctionBegin;
1192f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
11936a388c36SPeter Brune   linesearch->damping = damping;
11946a388c36SPeter Brune   PetscFunctionReturn(0);
11956a388c36SPeter Brune }
11966a388c36SPeter Brune 
119759405d5eSPeter Brune /*@
119859405d5eSPeter Brune    SNESLineSearchGetOrder - Gets the line search approximation order.
119959405d5eSPeter Brune 
120059405d5eSPeter Brune    Input Parameters:
120178bcb3b5SPeter Brune .  linesearch - linesearch context
120259405d5eSPeter Brune 
120359405d5eSPeter Brune    Output Parameters:
12048e557f58SPeter Brune .  order - The order
120559405d5eSPeter Brune 
120678bcb3b5SPeter Brune    Possible Values for order:
12073c7d6663SPeter Brune +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
12083c7d6663SPeter Brune .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
12093c7d6663SPeter Brune -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
121078bcb3b5SPeter Brune 
121159405d5eSPeter Brune    Level: intermediate
121259405d5eSPeter Brune 
1213db781477SPatrick Sanan .seealso: `SNESLineSearchSetOrder()`
121459405d5eSPeter Brune @*/
121559405d5eSPeter Brune 
12169371c9d4SSatish Balay PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch linesearch, PetscInt *order) {
121759405d5eSPeter Brune   PetscFunctionBegin;
121859405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1219534a8f05SLisandro Dalcin   PetscValidIntPointer(order, 2);
122059405d5eSPeter Brune   *order = linesearch->order;
122159405d5eSPeter Brune   PetscFunctionReturn(0);
122259405d5eSPeter Brune }
122359405d5eSPeter Brune 
122459405d5eSPeter Brune /*@
12251f8196a2SJed Brown    SNESLineSearchSetOrder - Sets the maximum order of the polynomial fit used in the line search
122659405d5eSPeter Brune 
122759405d5eSPeter Brune    Input Parameters:
1228a2b725a8SWilliam Gropp +  linesearch - linesearch context
1229a2b725a8SWilliam Gropp -  order - The damping parameter
123059405d5eSPeter Brune 
123159405d5eSPeter Brune    Level: intermediate
123259405d5eSPeter Brune 
123378bcb3b5SPeter Brune    Possible Values for order:
12343c7d6663SPeter Brune +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
12353c7d6663SPeter Brune .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
12363c7d6663SPeter Brune -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
123778bcb3b5SPeter Brune 
123859405d5eSPeter Brune    Notes:
123959405d5eSPeter Brune    Variable orders are supported by the following line searches:
124078bcb3b5SPeter Brune +  bt - cubic and quadratic
124178bcb3b5SPeter Brune -  cp - linear and quadratic
124259405d5eSPeter Brune 
1243db781477SPatrick Sanan .seealso: `SNESLineSearchGetOrder()`, `SNESLineSearchSetDamping()`
124459405d5eSPeter Brune @*/
12459371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch linesearch, PetscInt order) {
124659405d5eSPeter Brune   PetscFunctionBegin;
124759405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
124859405d5eSPeter Brune   linesearch->order = order;
124959405d5eSPeter Brune   PetscFunctionReturn(0);
125059405d5eSPeter Brune }
125159405d5eSPeter Brune 
1252f40b411bSPeter Brune /*@
1253f1c6b773SPeter Brune    SNESLineSearchGetNorms - Gets the norms for for X, Y, and F.
1254f40b411bSPeter Brune 
1255f899ff85SJose E. Roman    Input Parameter:
125678bcb3b5SPeter Brune .  linesearch - linesearch context
1257f40b411bSPeter Brune 
1258f40b411bSPeter Brune    Output Parameters:
1259f40b411bSPeter Brune +  xnorm - The norm of the current solution
1260f40b411bSPeter Brune .  fnorm - The norm of the current function
1261f40b411bSPeter Brune -  ynorm - The norm of the current update
1262f40b411bSPeter Brune 
1263cd7522eaSPeter Brune    Notes:
1264cd7522eaSPeter Brune    This function is mainly called from SNES implementations.
1265cd7522eaSPeter Brune 
126678bcb3b5SPeter Brune    Level: developer
1267f40b411bSPeter Brune 
1268db781477SPatrick Sanan .seealso: `SNESLineSearchSetNorms()` `SNESLineSearchGetVecs()`
1269f40b411bSPeter Brune @*/
12709371c9d4SSatish Balay PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal *xnorm, PetscReal *fnorm, PetscReal *ynorm) {
1271bf7f4e0aSPeter Brune   PetscFunctionBegin;
1272f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1273f5af7f23SKarl Rupp   if (xnorm) *xnorm = linesearch->xnorm;
1274f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
1275f5af7f23SKarl Rupp   if (ynorm) *ynorm = linesearch->ynorm;
1276bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1277bf7f4e0aSPeter Brune }
1278bf7f4e0aSPeter Brune 
1279f40b411bSPeter Brune /*@
1280f1c6b773SPeter Brune    SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F.
1281f40b411bSPeter Brune 
1282f40b411bSPeter Brune    Input Parameters:
128378bcb3b5SPeter Brune +  linesearch - linesearch context
1284f40b411bSPeter Brune .  xnorm - The norm of the current solution
1285f40b411bSPeter Brune .  fnorm - The norm of the current function
1286f40b411bSPeter Brune -  ynorm - The norm of the current update
1287f40b411bSPeter Brune 
128878bcb3b5SPeter Brune    Level: advanced
1289f40b411bSPeter Brune 
1290db781477SPatrick Sanan .seealso: `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()`
1291f40b411bSPeter Brune @*/
12929371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm) {
12936a388c36SPeter Brune   PetscFunctionBegin;
1294f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
12956a388c36SPeter Brune   linesearch->xnorm = xnorm;
12966a388c36SPeter Brune   linesearch->fnorm = fnorm;
12976a388c36SPeter Brune   linesearch->ynorm = ynorm;
12986a388c36SPeter Brune   PetscFunctionReturn(0);
12996a388c36SPeter Brune }
13006a388c36SPeter Brune 
1301f40b411bSPeter Brune /*@
1302f1c6b773SPeter Brune    SNESLineSearchComputeNorms - Computes the norms of X, F, and Y.
1303f40b411bSPeter Brune 
1304f40b411bSPeter Brune    Input Parameters:
130578bcb3b5SPeter Brune .  linesearch - linesearch context
1306f40b411bSPeter Brune 
1307f40b411bSPeter Brune    Options Database Keys:
13088e557f58SPeter Brune .   -snes_linesearch_norms - turn norm computation on or off
1309f40b411bSPeter Brune 
1310f40b411bSPeter Brune    Level: intermediate
1311f40b411bSPeter Brune 
1312db781477SPatrick Sanan .seealso: `SNESLineSearchGetNorms`, `SNESLineSearchSetNorms()`, `SNESLineSearchSetComputeNorms()`
1313f40b411bSPeter Brune @*/
13149371c9d4SSatish Balay PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch) {
13159bd66eb0SPeter Brune   SNES snes;
1316bf388a1fSBarry Smith 
13176a388c36SPeter Brune   PetscFunctionBegin;
13186a388c36SPeter Brune   if (linesearch->norms) {
13199bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
13209566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
13219566063dSJacob Faibussowitsch       PetscCall(VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
13229566063dSJacob Faibussowitsch       PetscCall(VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm));
13239566063dSJacob Faibussowitsch       PetscCall((*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm));
13249bd66eb0SPeter Brune     } else {
13259566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(linesearch->vec_func, NORM_2, &linesearch->fnorm));
13269566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
13279566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm));
13289566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(linesearch->vec_func, NORM_2, &linesearch->fnorm));
13299566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
13309566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(linesearch->vec_update, NORM_2, &linesearch->ynorm));
13316a388c36SPeter Brune     }
13329bd66eb0SPeter Brune   }
13336a388c36SPeter Brune   PetscFunctionReturn(0);
13346a388c36SPeter Brune }
13356a388c36SPeter Brune 
13366f263ca3SPeter Brune /*@
13376f263ca3SPeter Brune    SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search.
13386f263ca3SPeter Brune 
13396f263ca3SPeter Brune    Input Parameters:
134078bcb3b5SPeter Brune +  linesearch  - linesearch context
134178bcb3b5SPeter Brune -  flg  - indicates whether or not to compute norms
13426f263ca3SPeter Brune 
13436f263ca3SPeter Brune    Options Database Keys:
13440b00b554SBarry Smith .   -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic (none) linesearch
13456f263ca3SPeter Brune 
13466f263ca3SPeter Brune    Notes:
13470b00b554SBarry 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.
13486f263ca3SPeter Brune 
13496f263ca3SPeter Brune    Level: intermediate
13506f263ca3SPeter Brune 
1351db781477SPatrick Sanan .seealso: `SNESLineSearchGetNorms()`, `SNESLineSearchSetNorms()`, `SNESLineSearchComputeNorms()`, `SNESLINESEARCHBASIC`
13526f263ca3SPeter Brune @*/
13539371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg) {
13546f263ca3SPeter Brune   PetscFunctionBegin;
13556f263ca3SPeter Brune   linesearch->norms = flg;
13566f263ca3SPeter Brune   PetscFunctionReturn(0);
13576f263ca3SPeter Brune }
13586f263ca3SPeter Brune 
1359f40b411bSPeter Brune /*@
1360f1c6b773SPeter Brune    SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context
1361f40b411bSPeter Brune 
1362f899ff85SJose E. Roman    Input Parameter:
136378bcb3b5SPeter Brune .  linesearch - linesearch context
1364f40b411bSPeter Brune 
1365f40b411bSPeter Brune    Output Parameters:
13666232e825SPeter Brune +  X - Solution vector
13676232e825SPeter Brune .  F - Function vector
13686232e825SPeter Brune .  Y - Search direction vector
13696232e825SPeter Brune .  W - Solution work vector
13706232e825SPeter Brune -  G - Function work vector
13716232e825SPeter Brune 
13727bba9028SPeter Brune    Notes:
13737bba9028SPeter Brune    At the beginning of a line search application, X should contain a
13746232e825SPeter Brune    solution and the vector F the function computed at X.  At the end of the
13756232e825SPeter Brune    line search application, X should contain the new solution, and F the
13766232e825SPeter Brune    function evaluated at the new solution.
1377f40b411bSPeter Brune 
13782a7a6963SBarry Smith    These vectors are owned by the SNESLineSearch and should not be destroyed by the caller
13792a7a6963SBarry Smith 
138078bcb3b5SPeter Brune    Level: advanced
1381f40b411bSPeter Brune 
1382db781477SPatrick Sanan .seealso: `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()`
1383f40b411bSPeter Brune @*/
13849371c9d4SSatish Balay PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch, Vec *X, Vec *F, Vec *Y, Vec *W, Vec *G) {
13856a388c36SPeter Brune   PetscFunctionBegin;
1386f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
13876a388c36SPeter Brune   if (X) {
13886a388c36SPeter Brune     PetscValidPointer(X, 2);
13896a388c36SPeter Brune     *X = linesearch->vec_sol;
13906a388c36SPeter Brune   }
13916a388c36SPeter Brune   if (F) {
13926a388c36SPeter Brune     PetscValidPointer(F, 3);
13936a388c36SPeter Brune     *F = linesearch->vec_func;
13946a388c36SPeter Brune   }
13956a388c36SPeter Brune   if (Y) {
13966a388c36SPeter Brune     PetscValidPointer(Y, 4);
13976a388c36SPeter Brune     *Y = linesearch->vec_update;
13986a388c36SPeter Brune   }
13996a388c36SPeter Brune   if (W) {
14006a388c36SPeter Brune     PetscValidPointer(W, 5);
14016a388c36SPeter Brune     *W = linesearch->vec_sol_new;
14026a388c36SPeter Brune   }
14036a388c36SPeter Brune   if (G) {
14046a388c36SPeter Brune     PetscValidPointer(G, 6);
14056a388c36SPeter Brune     *G = linesearch->vec_func_new;
14066a388c36SPeter Brune   }
14076a388c36SPeter Brune   PetscFunctionReturn(0);
14086a388c36SPeter Brune }
14096a388c36SPeter Brune 
1410f40b411bSPeter Brune /*@
1411f1c6b773SPeter Brune    SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context
1412f40b411bSPeter Brune 
1413f40b411bSPeter Brune    Input Parameters:
141478bcb3b5SPeter Brune +  linesearch - linesearch context
14156232e825SPeter Brune .  X - Solution vector
14166232e825SPeter Brune .  F - Function vector
14176232e825SPeter Brune .  Y - Search direction vector
14186232e825SPeter Brune .  W - Solution work vector
14196232e825SPeter Brune -  G - Function work vector
1420f40b411bSPeter Brune 
142178bcb3b5SPeter Brune    Level: advanced
1422f40b411bSPeter Brune 
1423db781477SPatrick Sanan .seealso: `SNESLineSearchSetNorms()`, `SNESLineSearchGetVecs()`
1424f40b411bSPeter Brune @*/
14259371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch, Vec X, Vec F, Vec Y, Vec W, Vec G) {
14266a388c36SPeter Brune   PetscFunctionBegin;
1427f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
14286a388c36SPeter Brune   if (X) {
14296a388c36SPeter Brune     PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
14306a388c36SPeter Brune     linesearch->vec_sol = X;
14316a388c36SPeter Brune   }
14326a388c36SPeter Brune   if (F) {
14336a388c36SPeter Brune     PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
14346a388c36SPeter Brune     linesearch->vec_func = F;
14356a388c36SPeter Brune   }
14366a388c36SPeter Brune   if (Y) {
14376a388c36SPeter Brune     PetscValidHeaderSpecific(Y, VEC_CLASSID, 4);
14386a388c36SPeter Brune     linesearch->vec_update = Y;
14396a388c36SPeter Brune   }
14406a388c36SPeter Brune   if (W) {
14416a388c36SPeter Brune     PetscValidHeaderSpecific(W, VEC_CLASSID, 5);
14426a388c36SPeter Brune     linesearch->vec_sol_new = W;
14436a388c36SPeter Brune   }
14446a388c36SPeter Brune   if (G) {
14456a388c36SPeter Brune     PetscValidHeaderSpecific(G, VEC_CLASSID, 6);
14466a388c36SPeter Brune     linesearch->vec_func_new = G;
14476a388c36SPeter Brune   }
14486a388c36SPeter Brune   PetscFunctionReturn(0);
14496a388c36SPeter Brune }
14506a388c36SPeter Brune 
1451e7058c64SPeter Brune /*@C
1452f1c6b773SPeter Brune    SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1453e7058c64SPeter Brune    SNES options in the database.
1454e7058c64SPeter Brune 
1455cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
1456e7058c64SPeter Brune 
1457e7058c64SPeter Brune    Input Parameters:
1458e7058c64SPeter Brune +  snes - the SNES context
1459e7058c64SPeter Brune -  prefix - the prefix to prepend to all option names
1460e7058c64SPeter Brune 
1461e7058c64SPeter Brune    Notes:
1462e7058c64SPeter Brune    A hyphen (-) must NOT be given at the beginning of the prefix name.
1463e7058c64SPeter Brune    The first character of all runtime options is AUTOMATICALLY the hyphen.
1464e7058c64SPeter Brune 
1465e7058c64SPeter Brune    Level: advanced
1466e7058c64SPeter Brune 
1467db781477SPatrick Sanan .seealso: `SNESGetOptionsPrefix()`
1468e7058c64SPeter Brune @*/
14699371c9d4SSatish Balay PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch, const char prefix[]) {
1470e7058c64SPeter Brune   PetscFunctionBegin;
1471f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
14729566063dSJacob Faibussowitsch   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)linesearch, prefix));
1473e7058c64SPeter Brune   PetscFunctionReturn(0);
1474e7058c64SPeter Brune }
1475e7058c64SPeter Brune 
1476e7058c64SPeter Brune /*@C
1477f1c6b773SPeter Brune    SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all
1478f1c6b773SPeter Brune    SNESLineSearch options in the database.
1479e7058c64SPeter Brune 
1480e7058c64SPeter Brune    Not Collective
1481e7058c64SPeter Brune 
1482e7058c64SPeter Brune    Input Parameter:
1483cd7522eaSPeter Brune .  linesearch - the SNESLineSearch context
1484e7058c64SPeter Brune 
1485e7058c64SPeter Brune    Output Parameter:
1486e7058c64SPeter Brune .  prefix - pointer to the prefix string used
1487e7058c64SPeter Brune 
14888e557f58SPeter Brune    Notes:
14898e557f58SPeter Brune    On the fortran side, the user should pass in a string 'prefix' of
1490e7058c64SPeter Brune    sufficient length to hold the prefix.
1491e7058c64SPeter Brune 
1492e7058c64SPeter Brune    Level: advanced
1493e7058c64SPeter Brune 
1494db781477SPatrick Sanan .seealso: `SNESAppendOptionsPrefix()`
1495e7058c64SPeter Brune @*/
14969371c9d4SSatish Balay PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch, const char *prefix[]) {
1497e7058c64SPeter Brune   PetscFunctionBegin;
1498f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
14999566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)linesearch, prefix));
1500e7058c64SPeter Brune   PetscFunctionReturn(0);
1501e7058c64SPeter Brune }
1502bf7f4e0aSPeter Brune 
15038d359177SBarry Smith /*@C
1504fa0ddf94SBarry Smith    SNESLineSearchSetWorkVecs - Gets work vectors for the line search.
1505f40b411bSPeter Brune 
1506d8d19677SJose E. Roman    Input Parameters:
1507f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
1508f40b411bSPeter Brune -  nwork - the number of work vectors
1509f40b411bSPeter Brune 
1510f40b411bSPeter Brune    Level: developer
1511f40b411bSPeter Brune 
1512db781477SPatrick Sanan .seealso: `SNESSetWorkVecs()`
1513f40b411bSPeter Brune @*/
15149371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork) {
1515bf7f4e0aSPeter Brune   PetscFunctionBegin;
1516bf7f4e0aSPeter Brune   if (linesearch->vec_sol) {
15179566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work));
15188d359177SBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
1519bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1520bf7f4e0aSPeter Brune }
1521bf7f4e0aSPeter Brune 
1522f40b411bSPeter Brune /*@
1523422a814eSBarry Smith    SNESLineSearchGetReason - Gets the success/failure status of the last line search application
1524f40b411bSPeter Brune 
1525f40b411bSPeter Brune    Input Parameters:
152678bcb3b5SPeter Brune .  linesearch - linesearch context
1527f40b411bSPeter Brune 
1528f40b411bSPeter Brune    Output Parameters:
1529422a814eSBarry Smith .  result - The success or failure status
1530f40b411bSPeter Brune 
1531cd7522eaSPeter Brune    Notes:
1532c98378a5SBarry Smith    This is typically called after SNESLineSearchApply() in order to determine if the line-search failed
1533cd7522eaSPeter Brune    (and set the SNES convergence accordingly).
1534cd7522eaSPeter Brune 
1535f40b411bSPeter Brune    Level: intermediate
1536f40b411bSPeter Brune 
1537db781477SPatrick Sanan .seealso: `SNESLineSearchSetReason()`, `SNESLineSearchReason`
1538f40b411bSPeter Brune @*/
15399371c9d4SSatish Balay PetscErrorCode SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result) {
1540bf7f4e0aSPeter Brune   PetscFunctionBegin;
1541f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1542422a814eSBarry Smith   PetscValidPointer(result, 2);
1543422a814eSBarry Smith   *result = linesearch->result;
1544bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1545bf7f4e0aSPeter Brune }
1546bf7f4e0aSPeter Brune 
1547f40b411bSPeter Brune /*@
1548422a814eSBarry Smith    SNESLineSearchSetReason - Sets the success/failure status of the last line search application
1549f40b411bSPeter Brune 
1550f40b411bSPeter Brune    Input Parameters:
155178bcb3b5SPeter Brune +  linesearch - linesearch context
1552422a814eSBarry Smith -  result - The success or failure status
1553f40b411bSPeter Brune 
1554cd7522eaSPeter Brune    Notes:
1555cd7522eaSPeter Brune    This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set
1556cd7522eaSPeter Brune    the success or failure of the line search method.
1557cd7522eaSPeter Brune 
155878bcb3b5SPeter Brune    Level: developer
1559f40b411bSPeter Brune 
1560db781477SPatrick Sanan .seealso: `SNESLineSearchGetSResult()`
1561f40b411bSPeter Brune @*/
15629371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result) {
15636a388c36SPeter Brune   PetscFunctionBegin;
15645fd66863SKarl Rupp   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1565422a814eSBarry Smith   linesearch->result = result;
15666a388c36SPeter Brune   PetscFunctionReturn(0);
15676a388c36SPeter Brune }
15686a388c36SPeter Brune 
15699bd66eb0SPeter Brune /*@C
1570f1c6b773SPeter Brune    SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.
15719bd66eb0SPeter Brune 
15729bd66eb0SPeter Brune    Input Parameters:
15739bd66eb0SPeter Brune +  snes - nonlinear context obtained from SNESCreate()
15749bd66eb0SPeter Brune .  projectfunc - function for projecting the function to the bounds
15759bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
15769bd66eb0SPeter Brune 
15779bd66eb0SPeter Brune    Logically Collective on SNES
15789bd66eb0SPeter Brune 
15799bd66eb0SPeter Brune    Calling sequence of projectfunc:
15809bd66eb0SPeter Brune .vb
15819bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X)
15829bd66eb0SPeter Brune .ve
15839bd66eb0SPeter Brune 
15849bd66eb0SPeter Brune     Input parameters for projectfunc:
15859bd66eb0SPeter Brune +   snes - nonlinear context
15869bd66eb0SPeter Brune -   X - current solution
15879bd66eb0SPeter Brune 
1588cd7522eaSPeter Brune     Output parameters for projectfunc:
15899bd66eb0SPeter Brune .   X - Projected solution
15909bd66eb0SPeter Brune 
15919bd66eb0SPeter Brune    Calling sequence of normfunc:
15929bd66eb0SPeter Brune .vb
15939bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm)
15949bd66eb0SPeter Brune .ve
15959bd66eb0SPeter Brune 
1596cd7522eaSPeter Brune     Input parameters for normfunc:
15979bd66eb0SPeter Brune +   snes - nonlinear context
15989bd66eb0SPeter Brune .   X - current solution
15999bd66eb0SPeter Brune -   F - current residual
16009bd66eb0SPeter Brune 
1601cd7522eaSPeter Brune     Output parameters for normfunc:
16029bd66eb0SPeter Brune .   fnorm - VI-specific norm of the function
16039bd66eb0SPeter Brune 
1604cd7522eaSPeter Brune     Notes:
1605cd7522eaSPeter Brune     The VI solvers require projection of the solution to the feasible set.  projectfunc should implement this.
1606cd7522eaSPeter Brune 
1607cd7522eaSPeter Brune     The VI solvers require special evaluation of the function norm such that the norm is only calculated
1608cd7522eaSPeter Brune     on the inactive set.  This should be implemented by normfunc.
16099bd66eb0SPeter Brune 
16109bd66eb0SPeter Brune     Level: developer
16119bd66eb0SPeter Brune 
1612db781477SPatrick Sanan .seealso: `SNESLineSearchGetVIFunctions()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()`
16139bd66eb0SPeter Brune @*/
16149371c9d4SSatish Balay PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc) {
16159bd66eb0SPeter Brune   PetscFunctionBegin;
1616f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
16179bd66eb0SPeter Brune   if (projectfunc) linesearch->ops->viproject = projectfunc;
16189bd66eb0SPeter Brune   if (normfunc) linesearch->ops->vinorm = normfunc;
16199bd66eb0SPeter Brune   PetscFunctionReturn(0);
16209bd66eb0SPeter Brune }
16219bd66eb0SPeter Brune 
16229bd66eb0SPeter Brune /*@C
1623f1c6b773SPeter Brune    SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.
16249bd66eb0SPeter Brune 
1625f899ff85SJose E. Roman    Input Parameter:
1626907376e6SBarry Smith .  linesearch - the line search context, obtain with SNESGetLineSearch()
16279bd66eb0SPeter Brune 
16289bd66eb0SPeter Brune    Output Parameters:
16299bd66eb0SPeter Brune +  projectfunc - function for projecting the function to the bounds
16309bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
16319bd66eb0SPeter Brune 
16329bd66eb0SPeter Brune    Logically Collective on SNES
16339bd66eb0SPeter Brune 
16349bd66eb0SPeter Brune     Level: developer
16359bd66eb0SPeter Brune 
1636db781477SPatrick Sanan .seealso: `SNESLineSearchSetVIFunctions()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()`
16379bd66eb0SPeter Brune @*/
16389371c9d4SSatish Balay PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc) {
16399bd66eb0SPeter Brune   PetscFunctionBegin;
16409bd66eb0SPeter Brune   if (projectfunc) *projectfunc = linesearch->ops->viproject;
16419bd66eb0SPeter Brune   if (normfunc) *normfunc = linesearch->ops->vinorm;
16429bd66eb0SPeter Brune   PetscFunctionReturn(0);
16439bd66eb0SPeter Brune }
16449bd66eb0SPeter Brune 
1645bf7f4e0aSPeter Brune /*@C
16461c84c290SBarry Smith   SNESLineSearchRegister - See SNESLineSearchRegister()
1647bf7f4e0aSPeter Brune 
1648bf7f4e0aSPeter Brune   Level: advanced
1649bf7f4e0aSPeter Brune @*/
16509371c9d4SSatish Balay PetscErrorCode SNESLineSearchRegister(const char sname[], PetscErrorCode (*function)(SNESLineSearch)) {
1651bf7f4e0aSPeter Brune   PetscFunctionBegin;
16529566063dSJacob Faibussowitsch   PetscCall(SNESInitializePackage());
16539566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&SNESLineSearchList, sname, function));
1654bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1655bf7f4e0aSPeter Brune }
1656