xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision dbbe0bcd3f3a8fbab5a45420dc06f8387e5764c6)
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 @*/
32dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchMonitorCancel(SNESLineSearch ls)
33dcf2fd19SBarry Smith {
34dcf2fd19SBarry Smith   PetscInt       i;
35dcf2fd19SBarry Smith 
36dcf2fd19SBarry Smith   PetscFunctionBegin;
37dcf2fd19SBarry Smith   PetscValidHeaderSpecific(ls,SNESLINESEARCH_CLASSID,1);
38dcf2fd19SBarry Smith   for (i=0; i<ls->numbermonitors; i++) {
39dcf2fd19SBarry Smith     if (ls->monitordestroy[i]) {
409566063dSJacob Faibussowitsch       PetscCall((*ls->monitordestroy[i])(&ls->monitorcontext[i]));
41dcf2fd19SBarry Smith     }
42dcf2fd19SBarry Smith   }
43dcf2fd19SBarry Smith   ls->numbermonitors = 0;
44dcf2fd19SBarry Smith   PetscFunctionReturn(0);
45dcf2fd19SBarry Smith }
46dcf2fd19SBarry Smith 
47dcf2fd19SBarry Smith /*@
48dcf2fd19SBarry Smith    SNESLineSearchMonitor - runs the user provided monitor routines, if they exist
49dcf2fd19SBarry Smith 
50dcf2fd19SBarry Smith    Collective on SNES
51dcf2fd19SBarry Smith 
52dcf2fd19SBarry Smith    Input Parameters:
53dcf2fd19SBarry Smith .  ls - the linesearch object
54dcf2fd19SBarry Smith 
55dcf2fd19SBarry Smith    Notes:
56dcf2fd19SBarry Smith    This routine is called by the SNES implementations.
57dcf2fd19SBarry Smith    It does not typically need to be called by the user.
58dcf2fd19SBarry Smith 
59dcf2fd19SBarry Smith    Level: developer
60dcf2fd19SBarry Smith 
61db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchMonitorSet()`
62dcf2fd19SBarry Smith @*/
63dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchMonitor(SNESLineSearch ls)
64dcf2fd19SBarry Smith {
65dcf2fd19SBarry Smith   PetscInt       i,n = ls->numbermonitors;
66dcf2fd19SBarry Smith 
67dcf2fd19SBarry Smith   PetscFunctionBegin;
68dcf2fd19SBarry Smith   for (i=0; i<n; i++) {
699566063dSJacob Faibussowitsch     PetscCall((*ls->monitorftns[i])(ls,ls->monitorcontext[i]));
70dcf2fd19SBarry Smith   }
71dcf2fd19SBarry Smith   PetscFunctionReturn(0);
72dcf2fd19SBarry Smith }
73dcf2fd19SBarry Smith 
74dcf2fd19SBarry Smith /*@C
75dcf2fd19SBarry Smith    SNESLineSearchMonitorSet - Sets an ADDITIONAL function that is to be used at every
76dcf2fd19SBarry Smith    iteration of the nonlinear solver to display the iteration's
77dcf2fd19SBarry Smith    progress.
78dcf2fd19SBarry Smith 
79dcf2fd19SBarry Smith    Logically Collective on SNESLineSearch
80dcf2fd19SBarry Smith 
81dcf2fd19SBarry Smith    Input Parameters:
82dcf2fd19SBarry Smith +  ls - the SNESLineSearch context
83dcf2fd19SBarry Smith .  f - the monitor function
84dcf2fd19SBarry Smith .  mctx - [optional] user-defined context for private data for the
85dcf2fd19SBarry Smith           monitor routine (use NULL if no context is desired)
86dcf2fd19SBarry Smith -  monitordestroy - [optional] routine that frees monitor context
87dcf2fd19SBarry Smith           (may be NULL)
88dcf2fd19SBarry Smith 
89dcf2fd19SBarry Smith    Notes:
90dcf2fd19SBarry Smith    Several different monitoring routines may be set by calling
91dcf2fd19SBarry Smith    SNESLineSearchMonitorSet() multiple times; all will be called in the
92dcf2fd19SBarry Smith    order in which they were set.
93dcf2fd19SBarry Smith 
9495452b02SPatrick Sanan    Fortran Notes:
9595452b02SPatrick Sanan     Only a single monitor function can be set for each SNESLineSearch object
96dcf2fd19SBarry Smith 
97dcf2fd19SBarry Smith    Level: intermediate
98dcf2fd19SBarry Smith 
99db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchMonitorDefault()`, `SNESLineSearchMonitorCancel()`
100dcf2fd19SBarry Smith @*/
101dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchMonitorSet(SNESLineSearch ls,PetscErrorCode (*f)(SNESLineSearch,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
102dcf2fd19SBarry Smith {
10378064530SBarry Smith   PetscInt       i;
10478064530SBarry Smith   PetscBool      identical;
10578064530SBarry Smith 
106dcf2fd19SBarry Smith   PetscFunctionBegin;
107dcf2fd19SBarry Smith   PetscValidHeaderSpecific(ls,SNESLINESEARCH_CLASSID,1);
10878064530SBarry Smith   for (i=0; i<ls->numbermonitors;i++) {
1099566063dSJacob Faibussowitsch     PetscCall(PetscMonitorCompare((PetscErrorCode (*)(void))f,mctx,monitordestroy,(PetscErrorCode (*)(void))ls->monitorftns[i],ls->monitorcontext[i],ls->monitordestroy[i],&identical));
11078064530SBarry Smith     if (identical) PetscFunctionReturn(0);
11178064530SBarry Smith   }
1125f80ce2aSJacob Faibussowitsch   PetscCheck(ls->numbermonitors < MAXSNESLSMONITORS,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
113dcf2fd19SBarry Smith   ls->monitorftns[ls->numbermonitors]          = f;
114dcf2fd19SBarry Smith   ls->monitordestroy[ls->numbermonitors]   = monitordestroy;
115dcf2fd19SBarry Smith   ls->monitorcontext[ls->numbermonitors++] = (void*)mctx;
116dcf2fd19SBarry Smith   PetscFunctionReturn(0);
117dcf2fd19SBarry Smith }
118dcf2fd19SBarry Smith 
119dcf2fd19SBarry Smith /*@C
120dcf2fd19SBarry Smith    SNESLineSearchMonitorSolutionUpdate - Monitors each update a new function value the linesearch tries
121dcf2fd19SBarry Smith 
122dcf2fd19SBarry Smith    Collective on SNESLineSearch
123dcf2fd19SBarry Smith 
124dcf2fd19SBarry Smith    Input Parameters:
125dcf2fd19SBarry Smith +  ls - the SNES linesearch object
126d12e167eSBarry Smith -  vf - the context for the monitor, in this case it is an ASCII PetscViewer and format
127dcf2fd19SBarry Smith 
128dcf2fd19SBarry Smith    Level: intermediate
129dcf2fd19SBarry Smith 
130db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESMonitorSet()`, `SNESMonitorSolution()`
131dcf2fd19SBarry Smith @*/
132d12e167eSBarry Smith PetscErrorCode  SNESLineSearchMonitorSolutionUpdate(SNESLineSearch ls,PetscViewerAndFormat *vf)
133dcf2fd19SBarry Smith {
134d12e167eSBarry Smith   PetscViewer    viewer = vf->viewer;
135dcf2fd19SBarry Smith   Vec            Y,W,G;
136dcf2fd19SBarry Smith 
137dcf2fd19SBarry Smith   PetscFunctionBegin;
1389566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetVecs(ls,NULL,NULL,&Y,&W,&G));
1399566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(viewer,vf->format));
1409566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"LineSearch attempted update to solution \n"));
1419566063dSJacob Faibussowitsch   PetscCall(VecView(Y,viewer));
1429566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"LineSearch attempted new solution \n"));
1439566063dSJacob Faibussowitsch   PetscCall(VecView(W,viewer));
1449566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"LineSearch attempted updated function value\n"));
1459566063dSJacob Faibussowitsch   PetscCall(VecView(G,viewer));
1469566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(viewer));
147dcf2fd19SBarry Smith   PetscFunctionReturn(0);
148dcf2fd19SBarry Smith }
149dcf2fd19SBarry Smith 
150f40b411bSPeter Brune /*@
151cd7522eaSPeter Brune    SNESLineSearchCreate - Creates the line search context.
152f40b411bSPeter Brune 
153cd7522eaSPeter Brune    Logically Collective on Comm
154f40b411bSPeter Brune 
155f40b411bSPeter Brune    Input Parameters:
156cd7522eaSPeter Brune .  comm - MPI communicator for the line search (typically from the associated SNES context).
157f40b411bSPeter Brune 
158f40b411bSPeter Brune    Output Parameters:
1598e557f58SPeter Brune .  outlinesearch - the new linesearch context
160f40b411bSPeter Brune 
161162e0bf5SPeter Brune    Level: developer
162162e0bf5SPeter Brune 
163162e0bf5SPeter Brune    Notes:
164162e0bf5SPeter Brune    The preferred calling sequence for users is to use SNESGetLineSearch() to acquire the SNESLineSearch instance
165162e0bf5SPeter Brune    already associated with the SNES.  This function is for developer use.
166f40b411bSPeter Brune 
167db781477SPatrick Sanan .seealso: `LineSearchDestroy()`, `SNESGetLineSearch()`
168f40b411bSPeter Brune @*/
169f40b411bSPeter Brune 
170bf388a1fSBarry Smith PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch)
171bf388a1fSBarry Smith {
172f1c6b773SPeter Brune   SNESLineSearch linesearch;
173bf388a1fSBarry Smith 
174bf7f4e0aSPeter Brune   PetscFunctionBegin;
175ea5d4fccSPeter Brune   PetscValidPointer(outlinesearch,2);
1769566063dSJacob Faibussowitsch   PetscCall(SNESInitializePackage());
1770298fd71SBarry Smith   *outlinesearch = NULL;
178f5af7f23SKarl Rupp 
1799566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(linesearch,SNESLINESEARCH_CLASSID, "SNESLineSearch","Linesearch","SNESLineSearch",comm,SNESLineSearchDestroy,SNESLineSearchView));
180bf7f4e0aSPeter Brune 
1810298fd71SBarry Smith   linesearch->vec_sol_new  = NULL;
1820298fd71SBarry Smith   linesearch->vec_func_new = NULL;
1830298fd71SBarry Smith   linesearch->vec_sol      = NULL;
1840298fd71SBarry Smith   linesearch->vec_func     = NULL;
1850298fd71SBarry Smith   linesearch->vec_update   = NULL;
1869bd66eb0SPeter Brune 
187bf7f4e0aSPeter Brune   linesearch->lambda       = 1.0;
188bf7f4e0aSPeter Brune   linesearch->fnorm        = 1.0;
189bf7f4e0aSPeter Brune   linesearch->ynorm        = 1.0;
190bf7f4e0aSPeter Brune   linesearch->xnorm        = 1.0;
191422a814eSBarry Smith   linesearch->result       = SNES_LINESEARCH_SUCCEEDED;
192bf7f4e0aSPeter Brune   linesearch->norms        = PETSC_TRUE;
193bf7f4e0aSPeter Brune   linesearch->keeplambda   = PETSC_FALSE;
194bf7f4e0aSPeter Brune   linesearch->damping      = 1.0;
195bf7f4e0aSPeter Brune   linesearch->maxstep      = 1e8;
196bf7f4e0aSPeter Brune   linesearch->steptol      = 1e-12;
197516fe3c3SPeter Brune   linesearch->rtol         = 1e-8;
198516fe3c3SPeter Brune   linesearch->atol         = 1e-15;
199516fe3c3SPeter Brune   linesearch->ltol         = 1e-8;
2000298fd71SBarry Smith   linesearch->precheckctx  = NULL;
2010298fd71SBarry Smith   linesearch->postcheckctx = NULL;
20259405d5eSPeter Brune   linesearch->max_its      = 1;
203bf7f4e0aSPeter Brune   linesearch->setupcalled  = PETSC_FALSE;
2043add74b1SBarry Smith   linesearch->monitor      = NULL;
205bf7f4e0aSPeter Brune   *outlinesearch           = linesearch;
206bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
207bf7f4e0aSPeter Brune }
208bf7f4e0aSPeter Brune 
209f40b411bSPeter Brune /*@
21078bcb3b5SPeter Brune    SNESLineSearchSetUp - Prepares the line search for being applied by allocating
21178bcb3b5SPeter Brune    any required vectors.
212f40b411bSPeter Brune 
213cd7522eaSPeter Brune    Collective on SNESLineSearch
214f40b411bSPeter Brune 
215f40b411bSPeter Brune    Input Parameters:
216f40b411bSPeter Brune .  linesearch - The LineSearch instance.
217f40b411bSPeter Brune 
218cd7522eaSPeter Brune    Notes:
219f190f2fcSBarry Smith    For most cases, this needn't be called by users or outside of SNESLineSearchApply().
220cd7522eaSPeter Brune    The only current case where this is called outside of this is for the VI
22178bcb3b5SPeter Brune    solvers, which modify the solution and work vectors before the first call
222cd7522eaSPeter Brune    of SNESLineSearchApply, requiring the SNESLineSearch work vectors to be
223cd7522eaSPeter Brune    allocated upfront.
224cd7522eaSPeter Brune 
22578bcb3b5SPeter Brune    Level: advanced
226f40b411bSPeter Brune 
227db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchReset()`
228f40b411bSPeter Brune @*/
229f40b411bSPeter Brune 
230bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch)
231bf388a1fSBarry Smith {
232bf7f4e0aSPeter Brune   PetscFunctionBegin;
233bf7f4e0aSPeter Brune   if (!((PetscObject)linesearch)->type_name) {
2349566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC));
235bf7f4e0aSPeter Brune   }
236bf7f4e0aSPeter Brune   if (!linesearch->setupcalled) {
2379bd66eb0SPeter Brune     if (!linesearch->vec_sol_new) {
2389566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new));
2399bd66eb0SPeter Brune     }
2409bd66eb0SPeter Brune     if (!linesearch->vec_func_new) {
2419566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new));
2429bd66eb0SPeter Brune     }
243*dbbe0bcdSBarry Smith     if (linesearch->ops->setup) PetscUseTypeMethod(linesearch,setup);
2449566063dSJacob Faibussowitsch     if (!linesearch->ops->snesfunc) PetscCall(SNESLineSearchSetFunction(linesearch,SNESComputeFunction));
245bf7f4e0aSPeter Brune     linesearch->lambda      = linesearch->damping;
246bf7f4e0aSPeter Brune     linesearch->setupcalled = PETSC_TRUE;
247bf7f4e0aSPeter Brune   }
248bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
249bf7f4e0aSPeter Brune }
250bf7f4e0aSPeter Brune 
251f40b411bSPeter Brune /*@
252f190f2fcSBarry Smith    SNESLineSearchReset - Undoes the SNESLineSearchSetUp() and deletes any Vecs or Mats allocated by the line search.
253f40b411bSPeter Brune 
254f1c6b773SPeter Brune    Collective on SNESLineSearch
255f40b411bSPeter Brune 
256f40b411bSPeter Brune    Input Parameters:
257f40b411bSPeter Brune .  linesearch - The LineSearch instance.
258f40b411bSPeter Brune 
25995452b02SPatrick Sanan    Notes:
26095452b02SPatrick Sanan     Usually only called by SNESReset()
261f190f2fcSBarry Smith 
262f190f2fcSBarry Smith    Level: developer
263f40b411bSPeter Brune 
264db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchSetUp()`
265f40b411bSPeter Brune @*/
266f40b411bSPeter Brune 
267bf388a1fSBarry Smith PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch)
268bf388a1fSBarry Smith {
269bf7f4e0aSPeter Brune   PetscFunctionBegin;
270*dbbe0bcdSBarry Smith   if (linesearch->ops->reset) PetscUseTypeMethod(linesearch,reset);
271f5af7f23SKarl Rupp 
2729566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&linesearch->vec_sol_new));
2739566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&linesearch->vec_func_new));
274bf7f4e0aSPeter Brune 
2759566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(linesearch->nwork, &linesearch->work));
276f5af7f23SKarl Rupp 
277bf7f4e0aSPeter Brune   linesearch->nwork       = 0;
278bf7f4e0aSPeter Brune   linesearch->setupcalled = PETSC_FALSE;
279bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
280bf7f4e0aSPeter Brune }
281bf7f4e0aSPeter Brune 
282ed07d7d7SPeter Brune /*@C
283f190f2fcSBarry Smith    SNESLineSearchSetFunction - Sets the function evaluation used by the SNES line search
284ed07d7d7SPeter Brune 
285ed07d7d7SPeter Brune    Input Parameters:
286ed07d7d7SPeter Brune .  linesearch - the SNESLineSearch context
287f190f2fcSBarry Smith +  func       - function evaluation routine
288ed07d7d7SPeter Brune 
289ed07d7d7SPeter Brune    Level: developer
290ed07d7d7SPeter Brune 
29195452b02SPatrick Sanan    Notes:
29295452b02SPatrick Sanan     This is used internally by PETSc and not called by users
293f190f2fcSBarry Smith 
294db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESSetFunction()`
295ed07d7d7SPeter Brune @*/
296ed07d7d7SPeter Brune PetscErrorCode  SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES,Vec,Vec))
297ed07d7d7SPeter Brune {
298ed07d7d7SPeter Brune   PetscFunctionBegin;
299ed07d7d7SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
300ed07d7d7SPeter Brune   linesearch->ops->snesfunc = func;
301ed07d7d7SPeter Brune   PetscFunctionReturn(0);
302ed07d7d7SPeter Brune }
303ed07d7d7SPeter Brune 
30486d74e61SPeter Brune /*@C
305f190f2fcSBarry Smith    SNESLineSearchSetPreCheck - Sets a user function that is called after the initial search direction has been computed but
306df3898eeSBarry Smith          before the line search routine has been applied. Allows the user to adjust the result of (usually a linear solve) that
307f190f2fcSBarry Smith          determined the search direction.
30886d74e61SPeter Brune 
309f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
31086d74e61SPeter Brune 
31186d74e61SPeter Brune    Input Parameters:
312f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
31384238204SBarry Smith .  func - [optional] function evaluation routine, see SNESLineSearchPreCheck() for the calling sequence
314f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
31586d74e61SPeter Brune 
31686d74e61SPeter Brune    Level: intermediate
31786d74e61SPeter Brune 
318f0b84518SBarry Smith    Notes:
319f0b84518SBarry Smith    Use `SNESLineSearchSetPostCheck()` to change the step after the line search.
320f0b84518SBarry Smith    search is complete.
321f0b84518SBarry Smith 
322f0b84518SBarry Smith    Use `SNESVISetVariableBounds()` and `SNESVISetComputeVariableBounds()` to cause `SNES` to automatically control the ranges of variables allowed.
323f0b84518SBarry Smith 
324f0b84518SBarry Smith .seealso: `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()`,
325f0b84518SBarry Smith           `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()
326f0b84518SBarry Smith 
32786d74e61SPeter Brune @*/
328f190f2fcSBarry Smith PetscErrorCode  SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void *ctx)
32986d74e61SPeter Brune {
3309bd66eb0SPeter Brune   PetscFunctionBegin;
331f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
332f190f2fcSBarry Smith   if (func) linesearch->ops->precheck = func;
33386d74e61SPeter Brune   if (ctx) linesearch->precheckctx = ctx;
33486d74e61SPeter Brune   PetscFunctionReturn(0);
33586d74e61SPeter Brune }
33686d74e61SPeter Brune 
33786d74e61SPeter Brune /*@C
33853deb899SBarry Smith    SNESLineSearchGetPreCheck - Gets the pre-check function for the line search routine.
33986d74e61SPeter Brune 
340f899ff85SJose E. Roman    Input Parameter:
341f1c6b773SPeter Brune .  linesearch - the SNESLineSearch context
34286d74e61SPeter Brune 
34386d74e61SPeter Brune    Output Parameters:
34484238204SBarry Smith +  func       - [optional] function evaluation routine, see SNESLineSearchPreCheck() for calling sequence
345f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
34686d74e61SPeter Brune 
34786d74e61SPeter Brune    Level: intermediate
34886d74e61SPeter Brune 
349db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchSetPostCheck()`
35086d74e61SPeter Brune @*/
351f190f2fcSBarry Smith PetscErrorCode  SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void **ctx)
35286d74e61SPeter Brune {
3539bd66eb0SPeter Brune   PetscFunctionBegin;
354f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
355f190f2fcSBarry Smith   if (func) *func = linesearch->ops->precheck;
35686d74e61SPeter Brune   if (ctx) *ctx = linesearch->precheckctx;
35786d74e61SPeter Brune   PetscFunctionReturn(0);
35886d74e61SPeter Brune }
35986d74e61SPeter Brune 
36086d74e61SPeter Brune /*@C
361f190f2fcSBarry Smith    SNESLineSearchSetPostCheck - Sets a user function that is called after the line search has been applied to determine the step
362f190f2fcSBarry Smith        direction and length. Allows the user a chance to change or override the decision of the line search routine
36386d74e61SPeter Brune 
364f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
36586d74e61SPeter Brune 
36686d74e61SPeter Brune    Input Parameters:
367f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
36884238204SBarry Smith .  func - [optional] function evaluation routine, see SNESLineSearchPostCheck()  for the calling sequence
369f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
37086d74e61SPeter Brune 
37186d74e61SPeter Brune    Level: intermediate
37286d74e61SPeter Brune 
373f0b84518SBarry Smith    Notes:
374f0b84518SBarry Smith    Use `SNESLineSearchSetPreCheck()` to change the step before the line search.
375f0b84518SBarry Smith    search is complete.
376f0b84518SBarry Smith 
377f0b84518SBarry Smith    Use `SNESVISetVariableBounds()` and `SNESVISetComputeVariableBounds()` to cause `SNES` to automatically control the ranges of variables allowed.
378f0b84518SBarry Smith 
379f0b84518SBarry Smith .seealso: `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchGetPostCheck()`,
380f0b84518SBarry Smith           `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()
38186d74e61SPeter Brune @*/
382f190f2fcSBarry Smith PetscErrorCode  SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void *ctx)
38386d74e61SPeter Brune {
38486d74e61SPeter Brune   PetscFunctionBegin;
385f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
386f190f2fcSBarry Smith   if (func) linesearch->ops->postcheck = func;
38786d74e61SPeter Brune   if (ctx) linesearch->postcheckctx = ctx;
38886d74e61SPeter Brune   PetscFunctionReturn(0);
38986d74e61SPeter Brune }
39086d74e61SPeter Brune 
39186d74e61SPeter Brune /*@C
392f1c6b773SPeter Brune    SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine.
39386d74e61SPeter Brune 
394f899ff85SJose E. Roman    Input Parameter:
395f1c6b773SPeter Brune .  linesearch - the SNESLineSearch context
39686d74e61SPeter Brune 
39786d74e61SPeter Brune    Output Parameters:
39884238204SBarry Smith +  func - [optional] function evaluation routine, see for the calling sequence SNESLineSearchPostCheck()
399f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
40086d74e61SPeter Brune 
40186d74e61SPeter Brune    Level: intermediate
40286d74e61SPeter Brune 
403db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`
40486d74e61SPeter Brune @*/
405f190f2fcSBarry Smith PetscErrorCode  SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void **ctx)
40686d74e61SPeter Brune {
4079bd66eb0SPeter Brune   PetscFunctionBegin;
408f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
409f190f2fcSBarry Smith   if (func) *func = linesearch->ops->postcheck;
41086d74e61SPeter Brune   if (ctx) *ctx = linesearch->postcheckctx;
41186d74e61SPeter Brune   PetscFunctionReturn(0);
41286d74e61SPeter Brune }
41386d74e61SPeter Brune 
414f40b411bSPeter Brune /*@
415f1c6b773SPeter Brune    SNESLineSearchPreCheck - Prepares the line search for being applied.
416f40b411bSPeter Brune 
417cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
418f40b411bSPeter Brune 
419f40b411bSPeter Brune    Input Parameters:
4207b1df9c1SPeter Brune +  linesearch - The linesearch instance.
4217b1df9c1SPeter Brune .  X - The current solution
4227b1df9c1SPeter Brune -  Y - The step direction
423f40b411bSPeter Brune 
424f40b411bSPeter Brune    Output Parameters:
4258e557f58SPeter Brune .  changed - Indicator that the precheck routine has changed anything
426f40b411bSPeter Brune 
427f0b84518SBarry Smith    Level: advanced
428f40b411bSPeter Brune 
429f0b84518SBarry Smith .seealso: `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`,
430f0b84518SBarry Smith           `SNESLineSearchGetPostCheck()``
431f40b411bSPeter Brune @*/
4327b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed)
433bf7f4e0aSPeter Brune {
434bf7f4e0aSPeter Brune   PetscFunctionBegin;
435bf7f4e0aSPeter Brune   *changed = PETSC_FALSE;
4366b2b7091SBarry Smith   if (linesearch->ops->precheck) {
437*dbbe0bcdSBarry Smith     PetscUseTypeMethod(linesearch,precheck , X, Y, changed, linesearch->precheckctx);
43838bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed,4);
439bf7f4e0aSPeter Brune   }
440bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
441bf7f4e0aSPeter Brune }
442bf7f4e0aSPeter Brune 
443f40b411bSPeter Brune /*@
444ef46b1a6SStefano Zampini    SNESLineSearchPostCheck - Hook to modify step direction or updated solution after a successful linesearch
445f40b411bSPeter Brune 
446cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
447f40b411bSPeter Brune 
448f40b411bSPeter Brune    Input Parameters:
4497b1df9c1SPeter Brune +  linesearch - The linesearch context
4507b1df9c1SPeter Brune .  X - The last solution
4517b1df9c1SPeter Brune .  Y - The step direction
4527b1df9c1SPeter Brune -  W - The updated solution, W = X + lambda*Y for some lambda
453f40b411bSPeter Brune 
454f40b411bSPeter Brune    Output Parameters:
45578bcb3b5SPeter Brune +  changed_Y - Indicator if the direction Y has been changed.
45678bcb3b5SPeter Brune -  changed_W - Indicator if the new candidate solution W has been changed.
457f40b411bSPeter Brune 
458f190f2fcSBarry Smith    Level: developer
459f40b411bSPeter Brune 
460db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPrecheck()`, `SNESLineSearchGetPrecheck()`
461f40b411bSPeter Brune @*/
4627b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch,Vec X,Vec Y,Vec W,PetscBool *changed_Y,PetscBool *changed_W)
463bf7f4e0aSPeter Brune {
464bf7f4e0aSPeter Brune   PetscFunctionBegin;
465bf7f4e0aSPeter Brune   *changed_Y = PETSC_FALSE;
466bf7f4e0aSPeter Brune   *changed_W = PETSC_FALSE;
4676b2b7091SBarry Smith   if (linesearch->ops->postcheck) {
468*dbbe0bcdSBarry Smith     PetscUseTypeMethod(linesearch,postcheck ,X,Y,W,changed_Y,changed_W,linesearch->postcheckctx);
46938bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed_Y,5);
47038bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed_W,6);
47186d74e61SPeter Brune   }
47286d74e61SPeter Brune   PetscFunctionReturn(0);
47386d74e61SPeter Brune }
47486d74e61SPeter Brune 
47586d74e61SPeter Brune /*@C
47686d74e61SPeter Brune    SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration
47786d74e61SPeter Brune 
478cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
47986d74e61SPeter Brune 
4804165533cSJose E. Roman    Input Parameters:
48186d74e61SPeter Brune +  linesearch - linesearch context
48286d74e61SPeter Brune .  X - base state for this step
483907376e6SBarry Smith -  ctx - context for this function
48486d74e61SPeter Brune 
48597bb3fdcSJose E. Roman    Input/Output Parameter:
48697bb3fdcSJose E. Roman .  Y - correction, possibly modified
48797bb3fdcSJose E. Roman 
48897bb3fdcSJose E. Roman    Output Parameter:
48997bb3fdcSJose E. Roman .  changed - flag indicating that Y was modified
49086d74e61SPeter Brune 
49186d74e61SPeter Brune    Options Database Key:
492cd7522eaSPeter Brune +  -snes_linesearch_precheck_picard - activate this routine
493cd7522eaSPeter Brune -  -snes_linesearch_precheck_picard_angle - angle
49486d74e61SPeter Brune 
49586d74e61SPeter Brune    Level: advanced
49686d74e61SPeter Brune 
49786d74e61SPeter Brune    Notes:
49886d74e61SPeter Brune    This function should be passed to SNESLineSearchSetPreCheck()
49986d74e61SPeter Brune 
50086d74e61SPeter Brune    The justification for this method involves the linear convergence of a Picard iteration
50186d74e61SPeter Brune    so the Picard linearization should be provided in place of the "Jacobian". This correction
50286d74e61SPeter Brune    is generally not useful when using a Newton linearization.
50386d74e61SPeter Brune 
50486d74e61SPeter Brune    Reference:
50586d74e61SPeter Brune    Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology.
50686d74e61SPeter Brune 
507db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchSetPreCheck()`
50886d74e61SPeter Brune @*/
509f1c6b773SPeter Brune PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed,void *ctx)
51086d74e61SPeter Brune {
51186d74e61SPeter Brune   PetscReal      angle = *(PetscReal*)linesearch->precheckctx;
51286d74e61SPeter Brune   Vec            Ylast;
51386d74e61SPeter Brune   PetscScalar    dot;
51486d74e61SPeter Brune   PetscInt       iter;
51586d74e61SPeter Brune   PetscReal      ynorm,ylastnorm,theta,angle_radians;
51686d74e61SPeter Brune   SNES           snes;
51786d74e61SPeter Brune 
51886d74e61SPeter Brune   PetscFunctionBegin;
5199566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
5209566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject*)&Ylast));
52186d74e61SPeter Brune   if (!Ylast) {
5229566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(Y,&Ylast));
5239566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject)Ylast));
5249566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)Ylast));
52586d74e61SPeter Brune   }
5269566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes,&iter));
52786d74e61SPeter Brune   if (iter < 2) {
5289566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y,Ylast));
52986d74e61SPeter Brune     *changed = PETSC_FALSE;
53086d74e61SPeter Brune     PetscFunctionReturn(0);
53186d74e61SPeter Brune   }
53286d74e61SPeter Brune 
5339566063dSJacob Faibussowitsch   PetscCall(VecDot(Y,Ylast,&dot));
5349566063dSJacob Faibussowitsch   PetscCall(VecNorm(Y,NORM_2,&ynorm));
5359566063dSJacob Faibussowitsch   PetscCall(VecNorm(Ylast,NORM_2,&ylastnorm));
53686d74e61SPeter Brune   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
537255453a1SBarry Smith   theta         = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0));
53886d74e61SPeter Brune   angle_radians = angle * PETSC_PI / 180.;
53986d74e61SPeter Brune   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
54086d74e61SPeter Brune     /* Modify the step Y */
54186d74e61SPeter Brune     PetscReal alpha,ydiffnorm;
5429566063dSJacob Faibussowitsch     PetscCall(VecAXPY(Ylast,-1.0,Y));
5439566063dSJacob Faibussowitsch     PetscCall(VecNorm(Ylast,NORM_2,&ydiffnorm));
544f85e2ce2SBarry Smith     alpha = (ydiffnorm > .001*ylastnorm) ? ylastnorm / ydiffnorm : 1000.0;
5459566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y,Ylast));
5469566063dSJacob Faibussowitsch     PetscCall(VecScale(Y,alpha));
5479566063dSJacob 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));
548a47ec85fSBarry Smith     *changed = PETSC_TRUE;
54986d74e61SPeter Brune   } else {
5509566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes,"Angle %14.12e degrees exceeds threshold %14.12e, no correction applied\n",(double)(theta*180./PETSC_PI),(double)angle));
5519566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y,Ylast));
55286d74e61SPeter Brune     *changed = PETSC_FALSE;
553bf7f4e0aSPeter Brune   }
554bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
555bf7f4e0aSPeter Brune }
556bf7f4e0aSPeter Brune 
557f40b411bSPeter Brune /*@
558cd7522eaSPeter Brune    SNESLineSearchApply - Computes the line-search update.
559f40b411bSPeter Brune 
560f1c6b773SPeter Brune    Collective on SNESLineSearch
561f40b411bSPeter Brune 
562f40b411bSPeter Brune    Input Parameters:
5638e557f58SPeter Brune +  linesearch - The linesearch context
5648e557f58SPeter Brune -  Y - The search direction
565f40b411bSPeter Brune 
5666b867d5aSJose E. Roman    Input/Output Parameters:
5676b867d5aSJose E. Roman +  X - The current solution, on output the new solution
5686b867d5aSJose E. Roman .  F - The current function, on output the new function
5696b867d5aSJose E. Roman -  fnorm - The current norm, on output the new function norm
570f40b411bSPeter Brune 
571cd7522eaSPeter Brune    Options Database Keys:
5720b00b554SBarry Smith + -snes_linesearch_type - basic (or equivalently none), bt, l2, cp, nleqerr, shell
573dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches
5741fe24845SBarry Smith . -snes_linesearch_damping - The linesearch damping parameter, default is 1.0 (no damping)
5751fe24845SBarry Smith . -snes_linesearch_norms   - Turn on/off the linesearch norms computation (SNESLineSearchSetComputeNorms())
5763c7d6663SPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess
5773c7d6663SPeter Brune - -snes_linesearch_max_it - The number of iterations for iterative line searches
578cd7522eaSPeter Brune 
579cd7522eaSPeter Brune    Notes:
580cd7522eaSPeter Brune    This is typically called from within a SNESSolve() implementation in order to
581cd7522eaSPeter Brune    help with convergence of the nonlinear method.  Various SNES types use line searches
582cd7522eaSPeter Brune    in different ways, but the overarching theme is that a line search is used to determine
583cd7522eaSPeter Brune    an optimal damping parameter of a step at each iteration of the method.  Each
58484238204SBarry Smith    application of the line search may invoke SNESComputeFunction() several times, and
585cd7522eaSPeter Brune    therefore may be fairly expensive.
586cd7522eaSPeter Brune 
587f40b411bSPeter Brune    Level: Intermediate
588f40b411bSPeter Brune 
589db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchPreCheck()`, `SNESLineSearchPostCheck()`, `SNESSolve()`, `SNESComputeFunction()`, `SNESLineSearchSetComputeNorms()`,
590db781477SPatrick Sanan           `SNESLineSearchType`, `SNESLineSearchSetType()`
591f40b411bSPeter Brune @*/
592bf388a1fSBarry Smith PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y)
593bf388a1fSBarry Smith {
594bf388a1fSBarry Smith   PetscFunctionBegin;
595f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
596bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
597bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
598064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(Y,VEC_CLASSID,5);
599bf7f4e0aSPeter Brune 
600422a814eSBarry Smith   linesearch->result = SNES_LINESEARCH_SUCCEEDED;
601bf7f4e0aSPeter Brune 
602bf7f4e0aSPeter Brune   linesearch->vec_sol    = X;
603bf7f4e0aSPeter Brune   linesearch->vec_update = Y;
604bf7f4e0aSPeter Brune   linesearch->vec_func   = F;
605bf7f4e0aSPeter Brune 
6069566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchSetUp(linesearch));
607bf7f4e0aSPeter Brune 
608f5af7f23SKarl Rupp   if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */
609bf7f4e0aSPeter Brune 
610f5af7f23SKarl Rupp   if (fnorm) linesearch->fnorm = *fnorm;
611f5af7f23SKarl Rupp   else {
6129566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &linesearch->fnorm));
613bf7f4e0aSPeter Brune   }
614bf7f4e0aSPeter Brune 
6159566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(SNESLINESEARCH_Apply,linesearch,X,F,Y));
616bf7f4e0aSPeter Brune 
617*dbbe0bcdSBarry Smith   PetscUseTypeMethod(linesearch,apply);
618bf7f4e0aSPeter Brune 
6199566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(SNESLINESEARCH_Apply,linesearch,X,F,Y));
620bf7f4e0aSPeter Brune 
621f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
622bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
623bf7f4e0aSPeter Brune }
624bf7f4e0aSPeter Brune 
625f40b411bSPeter Brune /*@
626f1c6b773SPeter Brune    SNESLineSearchDestroy - Destroys the line search instance.
627f40b411bSPeter Brune 
628f1c6b773SPeter Brune    Collective on SNESLineSearch
629f40b411bSPeter Brune 
630f40b411bSPeter Brune    Input Parameters:
6318e557f58SPeter Brune .  linesearch - The linesearch context
632f40b411bSPeter Brune 
63384238204SBarry Smith    Level: developer
634f40b411bSPeter Brune 
635db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchReset()`, `SNESDestroy()`
636f40b411bSPeter Brune @*/
637bf388a1fSBarry Smith PetscErrorCode SNESLineSearchDestroy(SNESLineSearch * linesearch)
638bf388a1fSBarry Smith {
639bf7f4e0aSPeter Brune   PetscFunctionBegin;
640bf7f4e0aSPeter Brune   if (!*linesearch) PetscFunctionReturn(0);
641f1c6b773SPeter Brune   PetscValidHeaderSpecific((*linesearch),SNESLINESEARCH_CLASSID,1);
6429e5d0892SLisandro Dalcin   if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = NULL; PetscFunctionReturn(0);}
6439566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsViewOff((PetscObject)*linesearch));
6449566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchReset(*linesearch));
645f5af7f23SKarl Rupp   if ((*linesearch)->ops->destroy) (*linesearch)->ops->destroy(*linesearch);
6469566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&(*linesearch)->monitor));
6479566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchMonitorCancel((*linesearch)));
6489566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(linesearch));
649bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
650bf7f4e0aSPeter Brune }
651bf7f4e0aSPeter Brune 
652f40b411bSPeter Brune /*@
653dcf2fd19SBarry Smith    SNESLineSearchSetDefaultMonitor - Turns on/off printing useful information and debugging output about the line search.
654bf7f4e0aSPeter Brune 
655bf7f4e0aSPeter Brune    Input Parameters:
656dcf2fd19SBarry Smith +  linesearch - the linesearch object
657dcf2fd19SBarry Smith -  viewer - an ASCII PetscViewer or NULL to turn off monitor
658bf7f4e0aSPeter Brune 
659dcf2fd19SBarry Smith    Logically Collective on SNESLineSearch
660bf7f4e0aSPeter Brune 
661bf7f4e0aSPeter Brune    Options Database:
662dcf2fd19SBarry Smith .   -snes_linesearch_monitor [:filename] - enables the monitor
663bf7f4e0aSPeter Brune 
664bf7f4e0aSPeter Brune    Level: intermediate
665bf7f4e0aSPeter Brune 
666dcf2fd19SBarry Smith    Developer Note: This monitor is implemented differently than the other SNESLineSearchMonitors that are set with
667d12e167eSBarry Smith      SNESLineSearchMonitorSet() since it is called in many locations of the line search routines to display aspects of the
668d12e167eSBarry Smith      line search that are not visible to the other monitors.
669dcf2fd19SBarry Smith 
670db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchGetDefaultMonitor()`, `PetscViewer`, `SNESLineSearchSetMonitor()`
671bf7f4e0aSPeter Brune @*/
672dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchSetDefaultMonitor(SNESLineSearch linesearch, PetscViewer viewer)
673bf7f4e0aSPeter Brune {
674bf7f4e0aSPeter Brune   PetscFunctionBegin;
6759566063dSJacob Faibussowitsch   if (viewer) PetscCall(PetscObjectReference((PetscObject)viewer));
6769566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&linesearch->monitor));
677dcf2fd19SBarry Smith   linesearch->monitor = viewer;
678bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
679bf7f4e0aSPeter Brune }
680bf7f4e0aSPeter Brune 
681f40b411bSPeter Brune /*@
682dcf2fd19SBarry Smith    SNESLineSearchGetDefaultMonitor - Gets the PetscViewer instance for the line search monitor.
6836a388c36SPeter Brune 
684f190f2fcSBarry Smith    Input Parameter:
6858e557f58SPeter Brune .  linesearch - linesearch context
686f40b411bSPeter Brune 
687f190f2fcSBarry Smith    Output Parameter:
6888e557f58SPeter Brune .  monitor - monitor context
689f40b411bSPeter Brune 
690f40b411bSPeter Brune    Logically Collective on SNES
691f40b411bSPeter Brune 
692f40b411bSPeter Brune    Options Database Keys:
6938e557f58SPeter Brune .   -snes_linesearch_monitor - enables the monitor
694f40b411bSPeter Brune 
695f40b411bSPeter Brune    Level: intermediate
696f40b411bSPeter Brune 
697db781477SPatrick Sanan .seealso: `SNESGetLineSearch()`, `SNESLineSearchSetDefaultMonitor()`, `PetscViewer`
698f40b411bSPeter Brune @*/
699dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchGetDefaultMonitor(SNESLineSearch linesearch, PetscViewer *monitor)
7006a388c36SPeter Brune {
7016a388c36SPeter Brune   PetscFunctionBegin;
702f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
7036a388c36SPeter Brune   *monitor = linesearch->monitor;
7046a388c36SPeter Brune   PetscFunctionReturn(0);
7056a388c36SPeter Brune }
7066a388c36SPeter Brune 
707dcf2fd19SBarry Smith /*@C
708dcf2fd19SBarry Smith    SNESLineSearchMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
709dcf2fd19SBarry Smith 
710dcf2fd19SBarry Smith    Collective on SNESLineSearch
711dcf2fd19SBarry Smith 
712dcf2fd19SBarry Smith    Input Parameters:
713dcf2fd19SBarry Smith +  ls - LineSearch object you wish to monitor
714dcf2fd19SBarry Smith .  name - the monitor type one is seeking
715dcf2fd19SBarry Smith .  help - message indicating what monitoring is done
716dcf2fd19SBarry Smith .  manual - manual page for the monitor
717dcf2fd19SBarry Smith .  monitor - the monitor function
718dcf2fd19SBarry 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
719dcf2fd19SBarry Smith 
720dcf2fd19SBarry Smith    Level: developer
721dcf2fd19SBarry Smith 
722db781477SPatrick Sanan .seealso: `PetscOptionsGetViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
723db781477SPatrick Sanan           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
724db781477SPatrick Sanan           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
725db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
726c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
727db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
728db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
729dcf2fd19SBarry Smith @*/
730d12e167eSBarry Smith PetscErrorCode  SNESLineSearchMonitorSetFromOptions(SNESLineSearch ls,const char name[],const char help[], const char manual[],PetscErrorCode (*monitor)(SNESLineSearch,PetscViewerAndFormat*),PetscErrorCode (*monitorsetup)(SNESLineSearch,PetscViewerAndFormat*))
731dcf2fd19SBarry Smith {
732dcf2fd19SBarry Smith   PetscViewer       viewer;
733dcf2fd19SBarry Smith   PetscViewerFormat format;
734dcf2fd19SBarry Smith   PetscBool         flg;
735dcf2fd19SBarry Smith 
736dcf2fd19SBarry Smith   PetscFunctionBegin;
7379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)ls),((PetscObject) ls)->options,((PetscObject)ls)->prefix,name,&viewer,&format,&flg));
738dcf2fd19SBarry Smith   if (flg) {
739d12e167eSBarry Smith     PetscViewerAndFormat *vf;
7409566063dSJacob Faibussowitsch     PetscCall(PetscViewerAndFormatCreate(viewer,format,&vf));
7419566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)viewer));
7421baa6e33SBarry Smith     if (monitorsetup) PetscCall((*monitorsetup)(ls,vf));
7439566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchMonitorSet(ls,(PetscErrorCode (*)(SNESLineSearch,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy));
744dcf2fd19SBarry Smith   }
745dcf2fd19SBarry Smith   PetscFunctionReturn(0);
746dcf2fd19SBarry Smith }
747dcf2fd19SBarry Smith 
748f40b411bSPeter Brune /*@
749f1c6b773SPeter Brune    SNESLineSearchSetFromOptions - Sets options for the line search
750f40b411bSPeter Brune 
751f899ff85SJose E. Roman    Input Parameter:
7528e557f58SPeter Brune .  linesearch - linesearch context
753f40b411bSPeter Brune 
754f40b411bSPeter Brune    Options Database Keys:
7550b00b554SBarry Smith + -snes_linesearch_type <type> - basic (or equivalently none), bt, l2, cp, nleqerr, shell
7563c7d6663SPeter Brune . -snes_linesearch_order <order> - 1, 2, 3.  Most types only support certain orders (bt supports 2 or 3)
7571fe24845SBarry Smith . -snes_linesearch_norms   - Turn on/off the linesearch norms for the basic linesearch typem (SNESLineSearchSetComputeNorms())
75871eef1aeSPeter Brune . -snes_linesearch_minlambda - The minimum step length
7591a9b3a06SPeter Brune . -snes_linesearch_maxstep - The maximum step size
7601a9b3a06SPeter Brune . -snes_linesearch_rtol - Relative tolerance for iterative line searches
7611a9b3a06SPeter Brune . -snes_linesearch_atol - Absolute tolerance for iterative line searches
7621a9b3a06SPeter Brune . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches
7631a9b3a06SPeter Brune . -snes_linesearch_max_it - The number of iterations for iterative line searches
764dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches
765dcf2fd19SBarry Smith . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine
7668e557f58SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter
767cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess.
7681a9b3a06SPeter Brune . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method
769d8d34be6SBarry Smith - -snes_linesearch_precheck_picard_angle - Angle used in Picard precheck method
770f40b411bSPeter Brune 
771f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
772f40b411bSPeter Brune 
773f40b411bSPeter Brune    Level: intermediate
774f40b411bSPeter Brune 
775db781477SPatrick Sanan .seealso: `SNESLineSearchCreate()`, `SNESLineSearchSetOrder()`, `SNESLineSearchSetType()`, `SNESLineSearchSetTolerances()`, `SNESLineSearchSetDamping()`, `SNESLineSearchPreCheckPicard()`,
776db781477SPatrick Sanan           `SNESLineSearchType`, `SNESLineSearchSetComputeNorms()`
777f40b411bSPeter Brune @*/
778bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch)
779bf388a1fSBarry Smith {
7801a4f838cSPeter Brune   const char        *deft = SNESLINESEARCHBASIC;
781bf7f4e0aSPeter Brune   char              type[256];
782bf7f4e0aSPeter Brune   PetscBool         flg, set;
783dcf2fd19SBarry Smith   PetscViewer       viewer;
784bf388a1fSBarry Smith 
785bf7f4e0aSPeter Brune   PetscFunctionBegin;
7869566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchRegisterAll());
787bf7f4e0aSPeter Brune 
788d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)linesearch);
789f5af7f23SKarl Rupp   if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name;
7909566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-snes_linesearch_type","Linesearch type","SNESLineSearchSetType",SNESLineSearchList,deft,type,256,&flg));
791bf7f4e0aSPeter Brune   if (flg) {
7929566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetType(linesearch,type));
793bf7f4e0aSPeter Brune   } else if (!((PetscObject)linesearch)->type_name) {
7949566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetType(linesearch,deft));
795bf7f4e0aSPeter Brune   }
796bf7f4e0aSPeter Brune 
7979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)linesearch),((PetscObject) linesearch)->options,((PetscObject)linesearch)->prefix,"-snes_linesearch_monitor",&viewer,NULL,&set));
798dcf2fd19SBarry Smith   if (set) {
7999566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetDefaultMonitor(linesearch,viewer));
8009566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
801dcf2fd19SBarry Smith   }
8029566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchMonitorSetFromOptions(linesearch,"-snes_linesearch_monitor_solution_update","View correction at each iteration","SNESLineSearchMonitorSolutionUpdate",SNESLineSearchMonitorSolutionUpdate,NULL));
803bf7f4e0aSPeter Brune 
8041a9b3a06SPeter Brune   /* tolerances */
8059566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_minlambda","Minimum step length","SNESLineSearchSetTolerances",linesearch->steptol,&linesearch->steptol,NULL));
8069566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_maxstep","Maximum step size","SNESLineSearchSetTolerances",linesearch->maxstep,&linesearch->maxstep,NULL));
8079566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_rtol","Relative tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->rtol,&linesearch->rtol,NULL));
8089566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_atol","Absolute tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->atol,&linesearch->atol,NULL));
8099566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_ltol","Change in lambda tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->ltol,&linesearch->ltol,NULL));
8109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_linesearch_max_it","Maximum iterations for iterative line searches","SNESLineSearchSetTolerances",linesearch->max_its,&linesearch->max_its,NULL));
8117a35526eSPeter Brune 
8121a9b3a06SPeter Brune   /* damping parameters */
8139566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_damping","Line search damping and initial step guess","SNESLineSearchSetDamping",linesearch->damping,&linesearch->damping,NULL));
8141a9b3a06SPeter Brune 
8159566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_linesearch_keeplambda","Use previous lambda as damping","SNESLineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,NULL));
8161a9b3a06SPeter Brune 
8171a9b3a06SPeter Brune   /* precheck */
8189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_linesearch_precheck_picard","Use a correction that sometimes improves convergence of Picard iteration","SNESLineSearchPreCheckPicard",flg,&flg,&set));
8197a35526eSPeter Brune   if (set) {
8207a35526eSPeter Brune     if (flg) {
8217a35526eSPeter Brune       linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */
822f5af7f23SKarl Rupp 
823d0609cedSBarry 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));
8249566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&linesearch->precheck_picard_angle));
8257a35526eSPeter Brune     } else {
8269566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchSetPreCheck(linesearch,NULL,NULL));
8277a35526eSPeter Brune     }
8287a35526eSPeter Brune   }
8299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_linesearch_order","Order of approximation used in the line search","SNESLineSearchSetOrder",linesearch->order,&linesearch->order,NULL));
8309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_linesearch_norms","Compute final norms in line search","SNESLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,NULL));
8317a35526eSPeter Brune 
832*dbbe0bcdSBarry Smith   PetscTryTypeMethod(linesearch,setfromoptions,PetscOptionsObject);
8335a9b6599SPeter Brune 
834*dbbe0bcdSBarry Smith   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)linesearch,PetscOptionsObject));
835d0609cedSBarry Smith   PetscOptionsEnd();
836bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
837bf7f4e0aSPeter Brune }
838bf7f4e0aSPeter Brune 
839f40b411bSPeter Brune /*@
840f190f2fcSBarry Smith    SNESLineSearchView - Prints useful information about the line search
841f40b411bSPeter Brune 
842f40b411bSPeter Brune    Input Parameters:
8438e557f58SPeter Brune .  linesearch - linesearch context
844f40b411bSPeter Brune 
845f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
846f40b411bSPeter Brune 
847f40b411bSPeter Brune    Level: intermediate
848f40b411bSPeter Brune 
849db781477SPatrick Sanan .seealso: `SNESLineSearchCreate()`
850f40b411bSPeter Brune @*/
851bf388a1fSBarry Smith PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer)
852bf388a1fSBarry Smith {
8537f1410a3SPeter Brune   PetscBool      iascii;
854bf388a1fSBarry Smith 
855bf7f4e0aSPeter Brune   PetscFunctionBegin;
8567f1410a3SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
8577f1410a3SPeter Brune   if (!viewer) {
8589566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch),&viewer));
8597f1410a3SPeter Brune   }
8607f1410a3SPeter Brune   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
8617f1410a3SPeter Brune   PetscCheckSameComm(linesearch,1,viewer,2);
862f40b411bSPeter Brune 
8639566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
8647f1410a3SPeter Brune   if (iascii) {
8659566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)linesearch,viewer));
8669566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
867*dbbe0bcdSBarry Smith     PetscTryTypeMethod(linesearch,view ,viewer);
8689566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
8699566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep,(double)linesearch->steptol));
8709566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol,(double)linesearch->atol,(double)linesearch->ltol));
87163a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  maximum iterations=%" PetscInt_FMT "\n", linesearch->max_its));
8726b2b7091SBarry Smith     if (linesearch->ops->precheck) {
8736b2b7091SBarry Smith       if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) {
87463a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  using precheck step to speed up Picard convergence\n"));
8757f1410a3SPeter Brune       } else {
87663a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  using user-defined precheck step\n"));
8777f1410a3SPeter Brune       }
8787f1410a3SPeter Brune     }
8796b2b7091SBarry Smith     if (linesearch->ops->postcheck) {
88063a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  using user-defined postcheck step\n"));
8817f1410a3SPeter Brune     }
8827f1410a3SPeter Brune   }
883bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
884bf7f4e0aSPeter Brune }
885bf7f4e0aSPeter Brune 
886ea5d4fccSPeter Brune /*@C
887a80ff896SJed Brown    SNESLineSearchGetType - Gets the linesearch type
888a80ff896SJed Brown 
889a80ff896SJed Brown    Logically Collective on SNESLineSearch
890a80ff896SJed Brown 
891a80ff896SJed Brown    Input Parameters:
892a80ff896SJed Brown .  linesearch - linesearch context
893a80ff896SJed Brown 
894a80ff896SJed Brown    Output Parameters:
895a80ff896SJed Brown -  type - The type of line search, or NULL if not set
896a80ff896SJed Brown 
897a80ff896SJed Brown    Level: intermediate
898a80ff896SJed Brown 
899db781477SPatrick Sanan .seealso: `SNESLineSearchCreate()`, `SNESLineSearchType`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchSetType()`
900a80ff896SJed Brown @*/
901a80ff896SJed Brown PetscErrorCode SNESLineSearchGetType(SNESLineSearch linesearch, SNESLineSearchType *type)
902a80ff896SJed Brown {
903a80ff896SJed Brown   PetscFunctionBegin;
904a80ff896SJed Brown   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
905dadcf809SJacob Faibussowitsch   PetscValidPointer(type,2);
906a80ff896SJed Brown   *type = ((PetscObject)linesearch)->type_name;
907a80ff896SJed Brown   PetscFunctionReturn(0);
908a80ff896SJed Brown }
909a80ff896SJed Brown 
910a80ff896SJed Brown /*@C
911f1c6b773SPeter Brune    SNESLineSearchSetType - Sets the linesearch type
912f40b411bSPeter Brune 
913f190f2fcSBarry Smith    Logically Collective on SNESLineSearch
914f190f2fcSBarry Smith 
915f40b411bSPeter Brune    Input Parameters:
9168e557f58SPeter Brune +  linesearch - linesearch context
917f40b411bSPeter Brune -  type - The type of line search to be used
918f40b411bSPeter Brune 
919cd7522eaSPeter Brune    Available Types:
9200b00b554SBarry Smith +  SNESLINESEARCHBASIC - (or equivalently SNESLINESEARCHNONE) Simple damping line search, defaults to using the full Newton step
9211fe24845SBarry Smith .  SNESLINESEARCHBT - Backtracking line search over the L2 norm of the function
9221fe24845SBarry Smith .  SNESLINESEARCHL2 - Secant line search over the L2 norm of the function
9231fe24845SBarry Smith .  SNESLINESEARCHCP - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x)
9241fe24845SBarry Smith .  SNESLINESEARCHNLEQERR - Affine-covariant error-oriented linesearch
9251fe24845SBarry Smith -  SNESLINESEARCHSHELL - User provided SNESLineSearch implementation
9261fe24845SBarry Smith 
9271fe24845SBarry Smith    Options Database:
9280b00b554SBarry Smith .  -snes_linesearch_type <type> - basic (or equivalently none), bt, l2, cp, nleqerr, shell
929cd7522eaSPeter Brune 
930f40b411bSPeter Brune    Level: intermediate
931f40b411bSPeter Brune 
932db781477SPatrick Sanan .seealso: `SNESLineSearchCreate()`, `SNESLineSearchType`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchGetType()`
933f40b411bSPeter Brune @*/
93419fd82e9SBarry Smith PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type)
935bf7f4e0aSPeter Brune {
936bf7f4e0aSPeter Brune   PetscBool      match;
9375f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(SNESLineSearch);
938bf7f4e0aSPeter Brune 
939bf7f4e0aSPeter Brune   PetscFunctionBegin;
940f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
941bf7f4e0aSPeter Brune   PetscValidCharPointer(type,2);
942bf7f4e0aSPeter Brune 
9439566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)linesearch,type,&match));
944bf7f4e0aSPeter Brune   if (match) PetscFunctionReturn(0);
945bf7f4e0aSPeter Brune 
9469566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(SNESLineSearchList,type,&r));
9475f80ce2aSJacob Faibussowitsch   PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type);
948bf7f4e0aSPeter Brune   /* Destroy the previous private linesearch context */
949bf7f4e0aSPeter Brune   if (linesearch->ops->destroy) {
9509566063dSJacob Faibussowitsch     PetscCall((*(linesearch)->ops->destroy)(linesearch));
9510298fd71SBarry Smith     linesearch->ops->destroy = NULL;
952bf7f4e0aSPeter Brune   }
953f1c6b773SPeter Brune   /* Reinitialize function pointers in SNESLineSearchOps structure */
9549e5d0892SLisandro Dalcin   linesearch->ops->apply          = NULL;
9559e5d0892SLisandro Dalcin   linesearch->ops->view           = NULL;
9569e5d0892SLisandro Dalcin   linesearch->ops->setfromoptions = NULL;
9579e5d0892SLisandro Dalcin   linesearch->ops->destroy        = NULL;
958bf7f4e0aSPeter Brune 
9599566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)linesearch,type));
9609566063dSJacob Faibussowitsch   PetscCall((*r)(linesearch));
961bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
962bf7f4e0aSPeter Brune }
963bf7f4e0aSPeter Brune 
964f40b411bSPeter Brune /*@
96578bcb3b5SPeter Brune    SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation.
966f40b411bSPeter Brune 
967f40b411bSPeter Brune    Input Parameters:
9688e557f58SPeter Brune +  linesearch - linesearch context
969f40b411bSPeter Brune -  snes - The snes instance
970f40b411bSPeter Brune 
97178bcb3b5SPeter Brune    Level: developer
97278bcb3b5SPeter Brune 
97378bcb3b5SPeter Brune    Notes:
974f190f2fcSBarry Smith    This happens automatically when the line search is obtained/created with
9757601faf0SJed Brown    SNESGetLineSearch().  This routine is therefore mainly called within SNES
97678bcb3b5SPeter Brune    implementations.
977f40b411bSPeter Brune 
9788141a3b9SPeter Brune    Level: developer
979f40b411bSPeter Brune 
980db781477SPatrick Sanan .seealso: `SNESLineSearchGetSNES()`, `SNESLineSearchSetVecs()`, `SNES`
981f40b411bSPeter Brune @*/
982bf388a1fSBarry Smith PetscErrorCode  SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes)
983bf388a1fSBarry Smith {
984bf7f4e0aSPeter Brune   PetscFunctionBegin;
985f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
986bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,2);
987bf7f4e0aSPeter Brune   linesearch->snes = snes;
988bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
989bf7f4e0aSPeter Brune }
990bf7f4e0aSPeter Brune 
991f40b411bSPeter Brune /*@
9928141a3b9SPeter Brune    SNESLineSearchGetSNES - Gets the SNES instance associated with the line search.
9938141a3b9SPeter Brune    Having an associated SNES is necessary because most line search implementations must be able to
9948141a3b9SPeter Brune    evaluate the function using SNESComputeFunction() for the associated SNES.  This routine
9958141a3b9SPeter Brune    is used in the line search implementations when one must get this associated SNES instance.
996f40b411bSPeter Brune 
997f40b411bSPeter Brune    Input Parameters:
9988e557f58SPeter Brune .  linesearch - linesearch context
999f40b411bSPeter Brune 
1000f40b411bSPeter Brune    Output Parameters:
1001f40b411bSPeter Brune .  snes - The snes instance
1002f40b411bSPeter Brune 
10038141a3b9SPeter Brune    Level: developer
1004f40b411bSPeter Brune 
1005db781477SPatrick Sanan .seealso: `SNESLineSearchGetSNES()`, `SNESLineSearchSetVecs()`, `SNES`
1006f40b411bSPeter Brune @*/
1007bf388a1fSBarry Smith PetscErrorCode  SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes)
1008bf388a1fSBarry Smith {
1009bf7f4e0aSPeter Brune   PetscFunctionBegin;
1010f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
10116a388c36SPeter Brune   PetscValidPointer(snes,2);
1012bf7f4e0aSPeter Brune   *snes = linesearch->snes;
1013bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1014bf7f4e0aSPeter Brune }
1015bf7f4e0aSPeter Brune 
1016f40b411bSPeter Brune /*@
1017f1c6b773SPeter Brune    SNESLineSearchGetLambda - Gets the last linesearch steplength discovered.
1018f40b411bSPeter Brune 
1019f40b411bSPeter Brune    Input Parameters:
10208e557f58SPeter Brune .  linesearch - linesearch context
1021f40b411bSPeter Brune 
1022f40b411bSPeter Brune    Output Parameters:
1023cd7522eaSPeter Brune .  lambda - The last steplength computed during SNESLineSearchApply()
1024f40b411bSPeter Brune 
102578bcb3b5SPeter Brune    Level: advanced
102678bcb3b5SPeter Brune 
10278e557f58SPeter Brune    Notes:
10288e557f58SPeter Brune    This is useful in methods where the solver is ill-scaled and
102978bcb3b5SPeter Brune    requires some adaptive notion of the difference in scale between the
103078bcb3b5SPeter Brune    solution and the function.  For instance, SNESQN may be scaled by the
103178bcb3b5SPeter Brune    line search lambda using the argument -snes_qn_scaling ls.
103278bcb3b5SPeter Brune 
1033db781477SPatrick Sanan .seealso: `SNESLineSearchSetLambda()`, `SNESLineSearchGetDamping()`, `SNESLineSearchApply()`
1034f40b411bSPeter Brune @*/
1035f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda)
10366a388c36SPeter Brune {
10376a388c36SPeter Brune   PetscFunctionBegin;
1038f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1039534a8f05SLisandro Dalcin   PetscValidRealPointer(lambda, 2);
10406a388c36SPeter Brune   *lambda = linesearch->lambda;
10416a388c36SPeter Brune   PetscFunctionReturn(0);
10426a388c36SPeter Brune }
10436a388c36SPeter Brune 
1044f40b411bSPeter Brune /*@
1045f1c6b773SPeter Brune    SNESLineSearchSetLambda - Sets the linesearch steplength.
1046f40b411bSPeter Brune 
1047f40b411bSPeter Brune    Input Parameters:
10488e557f58SPeter Brune +  linesearch - linesearch context
1049f40b411bSPeter Brune -  lambda - The last steplength.
1050f40b411bSPeter Brune 
1051cd7522eaSPeter Brune    Notes:
1052f190f2fcSBarry Smith    This routine is typically used within implementations of SNESLineSearchApply()
1053cd7522eaSPeter Brune    to set the final steplength.  This routine (and SNESLineSearchGetLambda()) were
1054cd7522eaSPeter Brune    added in order to facilitate Quasi-Newton methods that use the previous steplength
1055cd7522eaSPeter Brune    as an inner scaling parameter.
1056cd7522eaSPeter Brune 
105778bcb3b5SPeter Brune    Level: advanced
1058f40b411bSPeter Brune 
1059db781477SPatrick Sanan .seealso: `SNESLineSearchGetLambda()`
1060f40b411bSPeter Brune @*/
1061f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
10626a388c36SPeter Brune {
10636a388c36SPeter Brune   PetscFunctionBegin;
1064f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
10656a388c36SPeter Brune   linesearch->lambda = lambda;
10666a388c36SPeter Brune   PetscFunctionReturn(0);
10676a388c36SPeter Brune }
10686a388c36SPeter Brune 
1069f40b411bSPeter Brune /*@
10703c7d6663SPeter Brune    SNESLineSearchGetTolerances - Gets the tolerances for the linesearch.  These include
107178bcb3b5SPeter Brune    tolerances for the relative and absolute change in the function norm, the change
107278bcb3b5SPeter Brune    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
107378bcb3b5SPeter Brune    and the maximum number of iterations the line search procedure may take.
1074f40b411bSPeter Brune 
1075f899ff85SJose E. Roman    Input Parameter:
10768e557f58SPeter Brune .  linesearch - linesearch context
1077f40b411bSPeter Brune 
1078f40b411bSPeter Brune    Output Parameters:
1079516fe3c3SPeter Brune +  steptol - The minimum steplength
10806cc8e53bSPeter Brune .  maxstep - The maximum steplength
1081516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
1082516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
1083516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
1084516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
1085f40b411bSPeter Brune 
108678bcb3b5SPeter Brune    Level: intermediate
108778bcb3b5SPeter Brune 
108878bcb3b5SPeter Brune    Notes:
108978bcb3b5SPeter Brune    Different line searches may implement these parameters slightly differently as
10903c7d6663SPeter Brune    the type requires.
1091516fe3c3SPeter Brune 
1092db781477SPatrick Sanan .seealso: `SNESLineSearchSetTolerances()`
1093f40b411bSPeter Brune @*/
1094f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
10956a388c36SPeter Brune {
10966a388c36SPeter Brune   PetscFunctionBegin;
1097f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1098516fe3c3SPeter Brune   if (steptol) {
1099534a8f05SLisandro Dalcin     PetscValidRealPointer(steptol, 2);
11006a388c36SPeter Brune     *steptol = linesearch->steptol;
1101516fe3c3SPeter Brune   }
1102516fe3c3SPeter Brune   if (maxstep) {
1103534a8f05SLisandro Dalcin     PetscValidRealPointer(maxstep, 3);
1104516fe3c3SPeter Brune     *maxstep = linesearch->maxstep;
1105516fe3c3SPeter Brune   }
1106516fe3c3SPeter Brune   if (rtol) {
1107534a8f05SLisandro Dalcin     PetscValidRealPointer(rtol, 4);
1108516fe3c3SPeter Brune     *rtol = linesearch->rtol;
1109516fe3c3SPeter Brune   }
1110516fe3c3SPeter Brune   if (atol) {
1111534a8f05SLisandro Dalcin     PetscValidRealPointer(atol, 5);
1112516fe3c3SPeter Brune     *atol = linesearch->atol;
1113516fe3c3SPeter Brune   }
1114516fe3c3SPeter Brune   if (ltol) {
1115534a8f05SLisandro Dalcin     PetscValidRealPointer(ltol, 6);
1116516fe3c3SPeter Brune     *ltol = linesearch->ltol;
1117516fe3c3SPeter Brune   }
1118516fe3c3SPeter Brune   if (max_its) {
1119534a8f05SLisandro Dalcin     PetscValidIntPointer(max_its, 7);
1120516fe3c3SPeter Brune     *max_its = linesearch->max_its;
1121516fe3c3SPeter Brune   }
11226a388c36SPeter Brune   PetscFunctionReturn(0);
11236a388c36SPeter Brune }
11246a388c36SPeter Brune 
1125f40b411bSPeter Brune /*@
11263c7d6663SPeter Brune    SNESLineSearchSetTolerances -  Gets the tolerances for the linesearch.  These include
112778bcb3b5SPeter Brune    tolerances for the relative and absolute change in the function norm, the change
112878bcb3b5SPeter Brune    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
112978bcb3b5SPeter Brune    and the maximum number of iterations the line search procedure may take.
1130f40b411bSPeter Brune 
1131f40b411bSPeter Brune    Input Parameters:
11328e557f58SPeter Brune +  linesearch - linesearch context
1133516fe3c3SPeter Brune .  steptol - The minimum steplength
11346cc8e53bSPeter Brune .  maxstep - The maximum steplength
1135516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
1136516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
1137516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
1138516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
1139f40b411bSPeter Brune 
114078bcb3b5SPeter Brune    Notes:
11413c7d6663SPeter Brune    The user may choose to not set any of the tolerances using PETSC_DEFAULT in
114278bcb3b5SPeter Brune    place of an argument.
1143f40b411bSPeter Brune 
114478bcb3b5SPeter Brune    Level: intermediate
1145516fe3c3SPeter Brune 
1146db781477SPatrick Sanan .seealso: `SNESLineSearchGetTolerances()`
1147f40b411bSPeter Brune @*/
1148f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its)
11496a388c36SPeter Brune {
11506a388c36SPeter Brune   PetscFunctionBegin;
1151f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1152d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,steptol,2);
1153d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,maxstep,3);
1154d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,rtol,4);
1155d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,atol,5);
1156d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,ltol,6);
1157d3952184SSatish Balay   PetscValidLogicalCollectiveInt(linesearch,max_its,7);
1158d3952184SSatish Balay 
1159d3952184SSatish Balay   if (steptol!= PETSC_DEFAULT) {
11605f80ce2aSJacob Faibussowitsch     PetscCheck(steptol >= 0.0,PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Minimum step length %14.12e must be non-negative",(double)steptol);
11616a388c36SPeter Brune     linesearch->steptol = steptol;
1162d3952184SSatish Balay   }
1163d3952184SSatish Balay 
1164d3952184SSatish Balay   if (maxstep!= PETSC_DEFAULT) {
11655f80ce2aSJacob Faibussowitsch     PetscCheck(maxstep >= 0.0,PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum step length %14.12e must be non-negative",(double)maxstep);
1166516fe3c3SPeter Brune     linesearch->maxstep = maxstep;
1167d3952184SSatish Balay   }
1168d3952184SSatish Balay 
1169d3952184SSatish Balay   if (rtol != PETSC_DEFAULT) {
11702061ca32SJunchao 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);
1171516fe3c3SPeter Brune     linesearch->rtol = rtol;
1172d3952184SSatish Balay   }
1173d3952184SSatish Balay 
1174d3952184SSatish Balay   if (atol != PETSC_DEFAULT) {
11755f80ce2aSJacob Faibussowitsch     PetscCheck(atol >= 0.0,PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %14.12e must be non-negative",(double)atol);
1176516fe3c3SPeter Brune     linesearch->atol = atol;
1177d3952184SSatish Balay   }
1178d3952184SSatish Balay 
1179d3952184SSatish Balay   if (ltol != PETSC_DEFAULT) {
11805f80ce2aSJacob Faibussowitsch     PetscCheck(ltol >= 0.0,PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Lambda tolerance %14.12e must be non-negative",(double)ltol);
1181516fe3c3SPeter Brune     linesearch->ltol = ltol;
1182d3952184SSatish Balay   }
1183d3952184SSatish Balay 
1184d3952184SSatish Balay   if (max_its != PETSC_DEFAULT) {
118563a3b9bcSJacob Faibussowitsch     PetscCheck(max_its >= 0,PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %" PetscInt_FMT " must be non-negative",max_its);
1186516fe3c3SPeter Brune     linesearch->max_its = max_its;
1187d3952184SSatish Balay   }
11886a388c36SPeter Brune   PetscFunctionReturn(0);
11896a388c36SPeter Brune }
11906a388c36SPeter Brune 
1191f40b411bSPeter Brune /*@
1192f1c6b773SPeter Brune    SNESLineSearchGetDamping - Gets the line search damping parameter.
1193f40b411bSPeter Brune 
1194f40b411bSPeter Brune    Input Parameters:
11958e557f58SPeter Brune .  linesearch - linesearch context
1196f40b411bSPeter Brune 
1197f40b411bSPeter Brune    Output Parameters:
11988e557f58SPeter Brune .  damping - The damping parameter
1199f40b411bSPeter Brune 
120078bcb3b5SPeter Brune    Level: advanced
1201f40b411bSPeter Brune 
1202db781477SPatrick Sanan .seealso: `SNESLineSearchGetStepTolerance()`, `SNESQN`
1203f40b411bSPeter Brune @*/
1204f40b411bSPeter Brune 
1205f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping)
12066a388c36SPeter Brune {
12076a388c36SPeter Brune   PetscFunctionBegin;
1208f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1209534a8f05SLisandro Dalcin   PetscValidRealPointer(damping, 2);
12106a388c36SPeter Brune   *damping = linesearch->damping;
12116a388c36SPeter Brune   PetscFunctionReturn(0);
12126a388c36SPeter Brune }
12136a388c36SPeter Brune 
1214f40b411bSPeter Brune /*@
1215fd292e60Sprj-    SNESLineSearchSetDamping - Sets the line search damping parameter.
1216f40b411bSPeter Brune 
1217f40b411bSPeter Brune    Input Parameters:
121803fd4120SBarry Smith +  linesearch - linesearch context
121903fd4120SBarry Smith -  damping - The damping parameter
1220f40b411bSPeter Brune 
122103fd4120SBarry Smith    Options Database:
122203fd4120SBarry Smith .   -snes_linesearch_damping
1223f40b411bSPeter Brune    Level: intermediate
1224f40b411bSPeter Brune 
1225cd7522eaSPeter Brune    Notes:
12260b00b554SBarry Smith    The basic (also known as the none) line search merely takes the update step scaled by the damping parameter.
1227cd7522eaSPeter Brune    The use of the damping parameter in the l2 and cp line searches is much more subtle;
122878bcb3b5SPeter Brune    it is used as a starting point in calculating the secant step. However, the eventual
1229cd7522eaSPeter Brune    step may be of greater length than the damping parameter.  In the bt line search it is
1230cd7522eaSPeter Brune    used as the maximum possible step length, as the bt line search only backtracks.
1231cd7522eaSPeter Brune 
1232db781477SPatrick Sanan .seealso: `SNESLineSearchGetDamping()`
1233f40b411bSPeter Brune @*/
1234f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping)
12356a388c36SPeter Brune {
12366a388c36SPeter Brune   PetscFunctionBegin;
1237f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
12386a388c36SPeter Brune   linesearch->damping = damping;
12396a388c36SPeter Brune   PetscFunctionReturn(0);
12406a388c36SPeter Brune }
12416a388c36SPeter Brune 
124259405d5eSPeter Brune /*@
124359405d5eSPeter Brune    SNESLineSearchGetOrder - Gets the line search approximation order.
124459405d5eSPeter Brune 
124559405d5eSPeter Brune    Input Parameters:
124678bcb3b5SPeter Brune .  linesearch - linesearch context
124759405d5eSPeter Brune 
124859405d5eSPeter Brune    Output Parameters:
12498e557f58SPeter Brune .  order - The order
125059405d5eSPeter Brune 
125178bcb3b5SPeter Brune    Possible Values for order:
12523c7d6663SPeter Brune +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
12533c7d6663SPeter Brune .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
12543c7d6663SPeter Brune -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
125578bcb3b5SPeter Brune 
125659405d5eSPeter Brune    Level: intermediate
125759405d5eSPeter Brune 
1258db781477SPatrick Sanan .seealso: `SNESLineSearchSetOrder()`
125959405d5eSPeter Brune @*/
126059405d5eSPeter Brune 
1261b000cd8dSPeter Brune PetscErrorCode  SNESLineSearchGetOrder(SNESLineSearch linesearch,PetscInt *order)
126259405d5eSPeter Brune {
126359405d5eSPeter Brune   PetscFunctionBegin;
126459405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1265534a8f05SLisandro Dalcin   PetscValidIntPointer(order, 2);
126659405d5eSPeter Brune   *order = linesearch->order;
126759405d5eSPeter Brune   PetscFunctionReturn(0);
126859405d5eSPeter Brune }
126959405d5eSPeter Brune 
127059405d5eSPeter Brune /*@
12711f8196a2SJed Brown    SNESLineSearchSetOrder - Sets the maximum order of the polynomial fit used in the line search
127259405d5eSPeter Brune 
127359405d5eSPeter Brune    Input Parameters:
1274a2b725a8SWilliam Gropp +  linesearch - linesearch context
1275a2b725a8SWilliam Gropp -  order - The damping parameter
127659405d5eSPeter Brune 
127759405d5eSPeter Brune    Level: intermediate
127859405d5eSPeter Brune 
127978bcb3b5SPeter Brune    Possible Values for order:
12803c7d6663SPeter Brune +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
12813c7d6663SPeter Brune .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
12823c7d6663SPeter Brune -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
128378bcb3b5SPeter Brune 
128459405d5eSPeter Brune    Notes:
128559405d5eSPeter Brune    Variable orders are supported by the following line searches:
128678bcb3b5SPeter Brune +  bt - cubic and quadratic
128778bcb3b5SPeter Brune -  cp - linear and quadratic
128859405d5eSPeter Brune 
1289db781477SPatrick Sanan .seealso: `SNESLineSearchGetOrder()`, `SNESLineSearchSetDamping()`
129059405d5eSPeter Brune @*/
1291b000cd8dSPeter Brune PetscErrorCode  SNESLineSearchSetOrder(SNESLineSearch linesearch,PetscInt order)
129259405d5eSPeter Brune {
129359405d5eSPeter Brune   PetscFunctionBegin;
129459405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
129559405d5eSPeter Brune   linesearch->order = order;
129659405d5eSPeter Brune   PetscFunctionReturn(0);
129759405d5eSPeter Brune }
129859405d5eSPeter Brune 
1299f40b411bSPeter Brune /*@
1300f1c6b773SPeter Brune    SNESLineSearchGetNorms - Gets the norms for for X, Y, and F.
1301f40b411bSPeter Brune 
1302f899ff85SJose E. Roman    Input Parameter:
130378bcb3b5SPeter Brune .  linesearch - linesearch context
1304f40b411bSPeter Brune 
1305f40b411bSPeter Brune    Output Parameters:
1306f40b411bSPeter Brune +  xnorm - The norm of the current solution
1307f40b411bSPeter Brune .  fnorm - The norm of the current function
1308f40b411bSPeter Brune -  ynorm - The norm of the current update
1309f40b411bSPeter Brune 
1310cd7522eaSPeter Brune    Notes:
1311cd7522eaSPeter Brune    This function is mainly called from SNES implementations.
1312cd7522eaSPeter Brune 
131378bcb3b5SPeter Brune    Level: developer
1314f40b411bSPeter Brune 
1315db781477SPatrick Sanan .seealso: `SNESLineSearchSetNorms()` `SNESLineSearchGetVecs()`
1316f40b411bSPeter Brune @*/
1317f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm)
1318bf7f4e0aSPeter Brune {
1319bf7f4e0aSPeter Brune   PetscFunctionBegin;
1320f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1321f5af7f23SKarl Rupp   if (xnorm) *xnorm = linesearch->xnorm;
1322f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
1323f5af7f23SKarl Rupp   if (ynorm) *ynorm = linesearch->ynorm;
1324bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1325bf7f4e0aSPeter Brune }
1326bf7f4e0aSPeter Brune 
1327f40b411bSPeter Brune /*@
1328f1c6b773SPeter Brune    SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F.
1329f40b411bSPeter Brune 
1330f40b411bSPeter Brune    Input Parameters:
133178bcb3b5SPeter Brune +  linesearch - linesearch context
1332f40b411bSPeter Brune .  xnorm - The norm of the current solution
1333f40b411bSPeter Brune .  fnorm - The norm of the current function
1334f40b411bSPeter Brune -  ynorm - The norm of the current update
1335f40b411bSPeter Brune 
133678bcb3b5SPeter Brune    Level: advanced
1337f40b411bSPeter Brune 
1338db781477SPatrick Sanan .seealso: `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()`
1339f40b411bSPeter Brune @*/
1340f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
13416a388c36SPeter Brune {
13426a388c36SPeter Brune   PetscFunctionBegin;
1343f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
13446a388c36SPeter Brune   linesearch->xnorm = xnorm;
13456a388c36SPeter Brune   linesearch->fnorm = fnorm;
13466a388c36SPeter Brune   linesearch->ynorm = ynorm;
13476a388c36SPeter Brune   PetscFunctionReturn(0);
13486a388c36SPeter Brune }
13496a388c36SPeter Brune 
1350f40b411bSPeter Brune /*@
1351f1c6b773SPeter Brune    SNESLineSearchComputeNorms - Computes the norms of X, F, and Y.
1352f40b411bSPeter Brune 
1353f40b411bSPeter Brune    Input Parameters:
135478bcb3b5SPeter Brune .  linesearch - linesearch context
1355f40b411bSPeter Brune 
1356f40b411bSPeter Brune    Options Database Keys:
13578e557f58SPeter Brune .   -snes_linesearch_norms - turn norm computation on or off
1358f40b411bSPeter Brune 
1359f40b411bSPeter Brune    Level: intermediate
1360f40b411bSPeter Brune 
1361db781477SPatrick Sanan .seealso: `SNESLineSearchGetNorms`, `SNESLineSearchSetNorms()`, `SNESLineSearchSetComputeNorms()`
1362f40b411bSPeter Brune @*/
1363f1c6b773SPeter Brune PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
13646a388c36SPeter Brune {
13659bd66eb0SPeter Brune   SNES           snes;
1366bf388a1fSBarry Smith 
13676a388c36SPeter Brune   PetscFunctionBegin;
13686a388c36SPeter Brune   if (linesearch->norms) {
13699bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
13709566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
13719566063dSJacob Faibussowitsch       PetscCall(VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
13729566063dSJacob Faibussowitsch       PetscCall(VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm));
13739566063dSJacob Faibussowitsch       PetscCall((*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm));
13749bd66eb0SPeter Brune     } else {
13759566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(linesearch->vec_func,   NORM_2, &linesearch->fnorm));
13769566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(linesearch->vec_sol,    NORM_2, &linesearch->xnorm));
13779566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm));
13789566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(linesearch->vec_func,     NORM_2, &linesearch->fnorm));
13799566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(linesearch->vec_sol,      NORM_2, &linesearch->xnorm));
13809566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(linesearch->vec_update,   NORM_2, &linesearch->ynorm));
13816a388c36SPeter Brune     }
13829bd66eb0SPeter Brune   }
13836a388c36SPeter Brune   PetscFunctionReturn(0);
13846a388c36SPeter Brune }
13856a388c36SPeter Brune 
13866f263ca3SPeter Brune /*@
13876f263ca3SPeter Brune    SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search.
13886f263ca3SPeter Brune 
13896f263ca3SPeter Brune    Input Parameters:
139078bcb3b5SPeter Brune +  linesearch  - linesearch context
139178bcb3b5SPeter Brune -  flg  - indicates whether or not to compute norms
13926f263ca3SPeter Brune 
13936f263ca3SPeter Brune    Options Database Keys:
13940b00b554SBarry Smith .   -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic (none) linesearch
13956f263ca3SPeter Brune 
13966f263ca3SPeter Brune    Notes:
13970b00b554SBarry 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.
13986f263ca3SPeter Brune 
13996f263ca3SPeter Brune    Level: intermediate
14006f263ca3SPeter Brune 
1401db781477SPatrick Sanan .seealso: `SNESLineSearchGetNorms()`, `SNESLineSearchSetNorms()`, `SNESLineSearchComputeNorms()`, `SNESLINESEARCHBASIC`
14026f263ca3SPeter Brune @*/
14036f263ca3SPeter Brune PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg)
14046f263ca3SPeter Brune {
14056f263ca3SPeter Brune   PetscFunctionBegin;
14066f263ca3SPeter Brune   linesearch->norms = flg;
14076f263ca3SPeter Brune   PetscFunctionReturn(0);
14086f263ca3SPeter Brune }
14096f263ca3SPeter Brune 
1410f40b411bSPeter Brune /*@
1411f1c6b773SPeter Brune    SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context
1412f40b411bSPeter Brune 
1413f899ff85SJose E. Roman    Input Parameter:
141478bcb3b5SPeter Brune .  linesearch - linesearch context
1415f40b411bSPeter Brune 
1416f40b411bSPeter Brune    Output Parameters:
14176232e825SPeter Brune +  X - Solution vector
14186232e825SPeter Brune .  F - Function vector
14196232e825SPeter Brune .  Y - Search direction vector
14206232e825SPeter Brune .  W - Solution work vector
14216232e825SPeter Brune -  G - Function work vector
14226232e825SPeter Brune 
14237bba9028SPeter Brune    Notes:
14247bba9028SPeter Brune    At the beginning of a line search application, X should contain a
14256232e825SPeter Brune    solution and the vector F the function computed at X.  At the end of the
14266232e825SPeter Brune    line search application, X should contain the new solution, and F the
14276232e825SPeter Brune    function evaluated at the new solution.
1428f40b411bSPeter Brune 
14292a7a6963SBarry Smith    These vectors are owned by the SNESLineSearch and should not be destroyed by the caller
14302a7a6963SBarry Smith 
143178bcb3b5SPeter Brune    Level: advanced
1432f40b411bSPeter Brune 
1433db781477SPatrick Sanan .seealso: `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()`
1434f40b411bSPeter Brune @*/
1435bf388a1fSBarry Smith PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G)
1436bf388a1fSBarry Smith {
14376a388c36SPeter Brune   PetscFunctionBegin;
1438f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
14396a388c36SPeter Brune   if (X) {
14406a388c36SPeter Brune     PetscValidPointer(X, 2);
14416a388c36SPeter Brune     *X = linesearch->vec_sol;
14426a388c36SPeter Brune   }
14436a388c36SPeter Brune   if (F) {
14446a388c36SPeter Brune     PetscValidPointer(F, 3);
14456a388c36SPeter Brune     *F = linesearch->vec_func;
14466a388c36SPeter Brune   }
14476a388c36SPeter Brune   if (Y) {
14486a388c36SPeter Brune     PetscValidPointer(Y, 4);
14496a388c36SPeter Brune     *Y = linesearch->vec_update;
14506a388c36SPeter Brune   }
14516a388c36SPeter Brune   if (W) {
14526a388c36SPeter Brune     PetscValidPointer(W, 5);
14536a388c36SPeter Brune     *W = linesearch->vec_sol_new;
14546a388c36SPeter Brune   }
14556a388c36SPeter Brune   if (G) {
14566a388c36SPeter Brune     PetscValidPointer(G, 6);
14576a388c36SPeter Brune     *G = linesearch->vec_func_new;
14586a388c36SPeter Brune   }
14596a388c36SPeter Brune   PetscFunctionReturn(0);
14606a388c36SPeter Brune }
14616a388c36SPeter Brune 
1462f40b411bSPeter Brune /*@
1463f1c6b773SPeter Brune    SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context
1464f40b411bSPeter Brune 
1465f40b411bSPeter Brune    Input Parameters:
146678bcb3b5SPeter Brune +  linesearch - linesearch context
14676232e825SPeter Brune .  X - Solution vector
14686232e825SPeter Brune .  F - Function vector
14696232e825SPeter Brune .  Y - Search direction vector
14706232e825SPeter Brune .  W - Solution work vector
14716232e825SPeter Brune -  G - Function work vector
1472f40b411bSPeter Brune 
147378bcb3b5SPeter Brune    Level: advanced
1474f40b411bSPeter Brune 
1475db781477SPatrick Sanan .seealso: `SNESLineSearchSetNorms()`, `SNESLineSearchGetVecs()`
1476f40b411bSPeter Brune @*/
1477bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G)
1478bf388a1fSBarry Smith {
14796a388c36SPeter Brune   PetscFunctionBegin;
1480f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
14816a388c36SPeter Brune   if (X) {
14826a388c36SPeter Brune     PetscValidHeaderSpecific(X,VEC_CLASSID,2);
14836a388c36SPeter Brune     linesearch->vec_sol = X;
14846a388c36SPeter Brune   }
14856a388c36SPeter Brune   if (F) {
14866a388c36SPeter Brune     PetscValidHeaderSpecific(F,VEC_CLASSID,3);
14876a388c36SPeter Brune     linesearch->vec_func = F;
14886a388c36SPeter Brune   }
14896a388c36SPeter Brune   if (Y) {
14906a388c36SPeter Brune     PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
14916a388c36SPeter Brune     linesearch->vec_update = Y;
14926a388c36SPeter Brune   }
14936a388c36SPeter Brune   if (W) {
14946a388c36SPeter Brune     PetscValidHeaderSpecific(W,VEC_CLASSID,5);
14956a388c36SPeter Brune     linesearch->vec_sol_new = W;
14966a388c36SPeter Brune   }
14976a388c36SPeter Brune   if (G) {
14986a388c36SPeter Brune     PetscValidHeaderSpecific(G,VEC_CLASSID,6);
14996a388c36SPeter Brune     linesearch->vec_func_new = G;
15006a388c36SPeter Brune   }
15016a388c36SPeter Brune   PetscFunctionReturn(0);
15026a388c36SPeter Brune }
15036a388c36SPeter Brune 
1504e7058c64SPeter Brune /*@C
1505f1c6b773SPeter Brune    SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1506e7058c64SPeter Brune    SNES options in the database.
1507e7058c64SPeter Brune 
1508cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
1509e7058c64SPeter Brune 
1510e7058c64SPeter Brune    Input Parameters:
1511e7058c64SPeter Brune +  snes - the SNES context
1512e7058c64SPeter Brune -  prefix - the prefix to prepend to all option names
1513e7058c64SPeter Brune 
1514e7058c64SPeter Brune    Notes:
1515e7058c64SPeter Brune    A hyphen (-) must NOT be given at the beginning of the prefix name.
1516e7058c64SPeter Brune    The first character of all runtime options is AUTOMATICALLY the hyphen.
1517e7058c64SPeter Brune 
1518e7058c64SPeter Brune    Level: advanced
1519e7058c64SPeter Brune 
1520db781477SPatrick Sanan .seealso: `SNESGetOptionsPrefix()`
1521e7058c64SPeter Brune @*/
1522f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[])
1523e7058c64SPeter Brune {
1524e7058c64SPeter Brune   PetscFunctionBegin;
1525f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
15269566063dSJacob Faibussowitsch   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix));
1527e7058c64SPeter Brune   PetscFunctionReturn(0);
1528e7058c64SPeter Brune }
1529e7058c64SPeter Brune 
1530e7058c64SPeter Brune /*@C
1531f1c6b773SPeter Brune    SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all
1532f1c6b773SPeter Brune    SNESLineSearch options in the database.
1533e7058c64SPeter Brune 
1534e7058c64SPeter Brune    Not Collective
1535e7058c64SPeter Brune 
1536e7058c64SPeter Brune    Input Parameter:
1537cd7522eaSPeter Brune .  linesearch - the SNESLineSearch context
1538e7058c64SPeter Brune 
1539e7058c64SPeter Brune    Output Parameter:
1540e7058c64SPeter Brune .  prefix - pointer to the prefix string used
1541e7058c64SPeter Brune 
15428e557f58SPeter Brune    Notes:
15438e557f58SPeter Brune    On the fortran side, the user should pass in a string 'prefix' of
1544e7058c64SPeter Brune    sufficient length to hold the prefix.
1545e7058c64SPeter Brune 
1546e7058c64SPeter Brune    Level: advanced
1547e7058c64SPeter Brune 
1548db781477SPatrick Sanan .seealso: `SNESAppendOptionsPrefix()`
1549e7058c64SPeter Brune @*/
1550f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[])
1551e7058c64SPeter Brune {
1552e7058c64SPeter Brune   PetscFunctionBegin;
1553f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
15549566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix));
1555e7058c64SPeter Brune   PetscFunctionReturn(0);
1556e7058c64SPeter Brune }
1557bf7f4e0aSPeter Brune 
15588d359177SBarry Smith /*@C
1559fa0ddf94SBarry Smith    SNESLineSearchSetWorkVecs - Gets work vectors for the line search.
1560f40b411bSPeter Brune 
1561d8d19677SJose E. Roman    Input Parameters:
1562f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
1563f40b411bSPeter Brune -  nwork - the number of work vectors
1564f40b411bSPeter Brune 
1565f40b411bSPeter Brune    Level: developer
1566f40b411bSPeter Brune 
1567db781477SPatrick Sanan .seealso: `SNESSetWorkVecs()`
1568f40b411bSPeter Brune @*/
1569fa0ddf94SBarry Smith PetscErrorCode  SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork)
1570bf7f4e0aSPeter Brune {
1571bf7f4e0aSPeter Brune   PetscFunctionBegin;
1572bf7f4e0aSPeter Brune   if (linesearch->vec_sol) {
15739566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work));
15748d359177SBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
1575bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1576bf7f4e0aSPeter Brune }
1577bf7f4e0aSPeter Brune 
1578f40b411bSPeter Brune /*@
1579422a814eSBarry Smith    SNESLineSearchGetReason - Gets the success/failure status of the last line search application
1580f40b411bSPeter Brune 
1581f40b411bSPeter Brune    Input Parameters:
158278bcb3b5SPeter Brune .  linesearch - linesearch context
1583f40b411bSPeter Brune 
1584f40b411bSPeter Brune    Output Parameters:
1585422a814eSBarry Smith .  result - The success or failure status
1586f40b411bSPeter Brune 
1587cd7522eaSPeter Brune    Notes:
1588c98378a5SBarry Smith    This is typically called after SNESLineSearchApply() in order to determine if the line-search failed
1589cd7522eaSPeter Brune    (and set the SNES convergence accordingly).
1590cd7522eaSPeter Brune 
1591f40b411bSPeter Brune    Level: intermediate
1592f40b411bSPeter Brune 
1593db781477SPatrick Sanan .seealso: `SNESLineSearchSetReason()`, `SNESLineSearchReason`
1594f40b411bSPeter Brune @*/
1595422a814eSBarry Smith PetscErrorCode  SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result)
1596bf7f4e0aSPeter Brune {
1597bf7f4e0aSPeter Brune   PetscFunctionBegin;
1598f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1599422a814eSBarry Smith   PetscValidPointer(result, 2);
1600422a814eSBarry Smith   *result = linesearch->result;
1601bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1602bf7f4e0aSPeter Brune }
1603bf7f4e0aSPeter Brune 
1604f40b411bSPeter Brune /*@
1605422a814eSBarry Smith    SNESLineSearchSetReason - Sets the success/failure status of the last line search application
1606f40b411bSPeter Brune 
1607f40b411bSPeter Brune    Input Parameters:
160878bcb3b5SPeter Brune +  linesearch - linesearch context
1609422a814eSBarry Smith -  result - The success or failure status
1610f40b411bSPeter Brune 
1611cd7522eaSPeter Brune    Notes:
1612cd7522eaSPeter Brune    This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set
1613cd7522eaSPeter Brune    the success or failure of the line search method.
1614cd7522eaSPeter Brune 
161578bcb3b5SPeter Brune    Level: developer
1616f40b411bSPeter Brune 
1617db781477SPatrick Sanan .seealso: `SNESLineSearchGetSResult()`
1618f40b411bSPeter Brune @*/
1619422a814eSBarry Smith PetscErrorCode  SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result)
16206a388c36SPeter Brune {
16216a388c36SPeter Brune   PetscFunctionBegin;
16225fd66863SKarl Rupp   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1623422a814eSBarry Smith   linesearch->result = result;
16246a388c36SPeter Brune   PetscFunctionReturn(0);
16256a388c36SPeter Brune }
16266a388c36SPeter Brune 
16279bd66eb0SPeter Brune /*@C
1628f1c6b773SPeter Brune    SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.
16299bd66eb0SPeter Brune 
16309bd66eb0SPeter Brune    Input Parameters:
16319bd66eb0SPeter Brune +  snes - nonlinear context obtained from SNESCreate()
16329bd66eb0SPeter Brune .  projectfunc - function for projecting the function to the bounds
16339bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
16349bd66eb0SPeter Brune 
16359bd66eb0SPeter Brune    Logically Collective on SNES
16369bd66eb0SPeter Brune 
16379bd66eb0SPeter Brune    Calling sequence of projectfunc:
16389bd66eb0SPeter Brune .vb
16399bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X)
16409bd66eb0SPeter Brune .ve
16419bd66eb0SPeter Brune 
16429bd66eb0SPeter Brune     Input parameters for projectfunc:
16439bd66eb0SPeter Brune +   snes - nonlinear context
16449bd66eb0SPeter Brune -   X - current solution
16459bd66eb0SPeter Brune 
1646cd7522eaSPeter Brune     Output parameters for projectfunc:
16479bd66eb0SPeter Brune .   X - Projected solution
16489bd66eb0SPeter Brune 
16499bd66eb0SPeter Brune    Calling sequence of normfunc:
16509bd66eb0SPeter Brune .vb
16519bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm)
16529bd66eb0SPeter Brune .ve
16539bd66eb0SPeter Brune 
1654cd7522eaSPeter Brune     Input parameters for normfunc:
16559bd66eb0SPeter Brune +   snes - nonlinear context
16569bd66eb0SPeter Brune .   X - current solution
16579bd66eb0SPeter Brune -   F - current residual
16589bd66eb0SPeter Brune 
1659cd7522eaSPeter Brune     Output parameters for normfunc:
16609bd66eb0SPeter Brune .   fnorm - VI-specific norm of the function
16619bd66eb0SPeter Brune 
1662cd7522eaSPeter Brune     Notes:
1663cd7522eaSPeter Brune     The VI solvers require projection of the solution to the feasible set.  projectfunc should implement this.
1664cd7522eaSPeter Brune 
1665cd7522eaSPeter Brune     The VI solvers require special evaluation of the function norm such that the norm is only calculated
1666cd7522eaSPeter Brune     on the inactive set.  This should be implemented by normfunc.
16679bd66eb0SPeter Brune 
16689bd66eb0SPeter Brune     Level: developer
16699bd66eb0SPeter Brune 
1670db781477SPatrick Sanan .seealso: `SNESLineSearchGetVIFunctions()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()`
16719bd66eb0SPeter Brune @*/
167225acbd8eSLisandro Dalcin PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc)
16739bd66eb0SPeter Brune {
16749bd66eb0SPeter Brune   PetscFunctionBegin;
1675f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
16769bd66eb0SPeter Brune   if (projectfunc) linesearch->ops->viproject = projectfunc;
16779bd66eb0SPeter Brune   if (normfunc) linesearch->ops->vinorm = normfunc;
16789bd66eb0SPeter Brune   PetscFunctionReturn(0);
16799bd66eb0SPeter Brune }
16809bd66eb0SPeter Brune 
16819bd66eb0SPeter Brune /*@C
1682f1c6b773SPeter Brune    SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.
16839bd66eb0SPeter Brune 
1684f899ff85SJose E. Roman    Input Parameter:
1685907376e6SBarry Smith .  linesearch - the line search context, obtain with SNESGetLineSearch()
16869bd66eb0SPeter Brune 
16879bd66eb0SPeter Brune    Output Parameters:
16889bd66eb0SPeter Brune +  projectfunc - function for projecting the function to the bounds
16899bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
16909bd66eb0SPeter Brune 
16919bd66eb0SPeter Brune    Logically Collective on SNES
16929bd66eb0SPeter Brune 
16939bd66eb0SPeter Brune     Level: developer
16949bd66eb0SPeter Brune 
1695db781477SPatrick Sanan .seealso: `SNESLineSearchSetVIFunctions()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()`
16969bd66eb0SPeter Brune @*/
169725acbd8eSLisandro Dalcin PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc)
16989bd66eb0SPeter Brune {
16999bd66eb0SPeter Brune   PetscFunctionBegin;
17009bd66eb0SPeter Brune   if (projectfunc) *projectfunc = linesearch->ops->viproject;
17019bd66eb0SPeter Brune   if (normfunc) *normfunc = linesearch->ops->vinorm;
17029bd66eb0SPeter Brune   PetscFunctionReturn(0);
17039bd66eb0SPeter Brune }
17049bd66eb0SPeter Brune 
1705bf7f4e0aSPeter Brune /*@C
17061c84c290SBarry Smith   SNESLineSearchRegister - See SNESLineSearchRegister()
1707bf7f4e0aSPeter Brune 
1708bf7f4e0aSPeter Brune   Level: advanced
1709bf7f4e0aSPeter Brune @*/
1710bdf89e91SBarry Smith PetscErrorCode  SNESLineSearchRegister(const char sname[],PetscErrorCode (*function)(SNESLineSearch))
1711bf7f4e0aSPeter Brune {
1712bf7f4e0aSPeter Brune   PetscFunctionBegin;
17139566063dSJacob Faibussowitsch   PetscCall(SNESInitializePackage());
17149566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&SNESLineSearchList,sname,function));
1715bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1716bf7f4e0aSPeter Brune }
1717