xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision ef46b1a67e276116c83b5d4ce8efc2932ea4fc0a)
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 
3084238204SBarry Smith .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 
6184238204SBarry Smith .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 
9984238204SBarry Smith .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 
13084238204SBarry Smith .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 
167162e0bf5SPeter Brune .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 
22784238204SBarry Smith .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     }
243bf7f4e0aSPeter Brune     if (linesearch->ops->setup) {
2449566063dSJacob Faibussowitsch       PetscCall((*linesearch->ops->setup)(linesearch));
245bf7f4e0aSPeter Brune     }
2469566063dSJacob Faibussowitsch     if (!linesearch->ops->snesfunc) PetscCall(SNESLineSearchSetFunction(linesearch,SNESComputeFunction));
247bf7f4e0aSPeter Brune     linesearch->lambda      = linesearch->damping;
248bf7f4e0aSPeter Brune     linesearch->setupcalled = PETSC_TRUE;
249bf7f4e0aSPeter Brune   }
250bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
251bf7f4e0aSPeter Brune }
252bf7f4e0aSPeter Brune 
253f40b411bSPeter Brune /*@
254f190f2fcSBarry Smith    SNESLineSearchReset - Undoes the SNESLineSearchSetUp() and deletes any Vecs or Mats allocated by the line search.
255f40b411bSPeter Brune 
256f1c6b773SPeter Brune    Collective on SNESLineSearch
257f40b411bSPeter Brune 
258f40b411bSPeter Brune    Input Parameters:
259f40b411bSPeter Brune .  linesearch - The LineSearch instance.
260f40b411bSPeter Brune 
26195452b02SPatrick Sanan    Notes:
26295452b02SPatrick Sanan     Usually only called by SNESReset()
263f190f2fcSBarry Smith 
264f190f2fcSBarry Smith    Level: developer
265f40b411bSPeter Brune 
26684238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchSetUp()
267f40b411bSPeter Brune @*/
268f40b411bSPeter Brune 
269bf388a1fSBarry Smith PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch)
270bf388a1fSBarry Smith {
271bf7f4e0aSPeter Brune   PetscFunctionBegin;
2729566063dSJacob Faibussowitsch   if (linesearch->ops->reset) PetscCall((*linesearch->ops->reset)(linesearch));
273f5af7f23SKarl Rupp 
2749566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&linesearch->vec_sol_new));
2759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&linesearch->vec_func_new));
276bf7f4e0aSPeter Brune 
2779566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(linesearch->nwork, &linesearch->work));
278f5af7f23SKarl Rupp 
279bf7f4e0aSPeter Brune   linesearch->nwork       = 0;
280bf7f4e0aSPeter Brune   linesearch->setupcalled = PETSC_FALSE;
281bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
282bf7f4e0aSPeter Brune }
283bf7f4e0aSPeter Brune 
284ed07d7d7SPeter Brune /*@C
285f190f2fcSBarry Smith    SNESLineSearchSetFunction - Sets the function evaluation used by the SNES line search
286ed07d7d7SPeter Brune 
287ed07d7d7SPeter Brune    Input Parameters:
288ed07d7d7SPeter Brune .  linesearch - the SNESLineSearch context
289f190f2fcSBarry Smith +  func       - function evaluation routine
290ed07d7d7SPeter Brune 
291ed07d7d7SPeter Brune    Level: developer
292ed07d7d7SPeter Brune 
29395452b02SPatrick Sanan    Notes:
29495452b02SPatrick Sanan     This is used internally by PETSc and not called by users
295f190f2fcSBarry Smith 
29684238204SBarry Smith .seealso: SNESGetLineSearch(), SNESSetFunction()
297ed07d7d7SPeter Brune @*/
298ed07d7d7SPeter Brune PetscErrorCode  SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES,Vec,Vec))
299ed07d7d7SPeter Brune {
300ed07d7d7SPeter Brune   PetscFunctionBegin;
301ed07d7d7SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
302ed07d7d7SPeter Brune   linesearch->ops->snesfunc = func;
303ed07d7d7SPeter Brune   PetscFunctionReturn(0);
304ed07d7d7SPeter Brune }
305ed07d7d7SPeter Brune 
30686d74e61SPeter Brune /*@C
307f190f2fcSBarry Smith    SNESLineSearchSetPreCheck - Sets a user function that is called after the initial search direction has been computed but
308df3898eeSBarry Smith          before the line search routine has been applied. Allows the user to adjust the result of (usually a linear solve) that
309f190f2fcSBarry Smith          determined the search direction.
31086d74e61SPeter Brune 
311f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
31286d74e61SPeter Brune 
31386d74e61SPeter Brune    Input Parameters:
314f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
31584238204SBarry Smith .  func - [optional] function evaluation routine, see SNESLineSearchPreCheck() for the calling sequence
316f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
31786d74e61SPeter Brune 
31886d74e61SPeter Brune    Level: intermediate
31986d74e61SPeter Brune 
32084238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
32186d74e61SPeter Brune @*/
322f190f2fcSBarry Smith PetscErrorCode  SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void *ctx)
32386d74e61SPeter Brune {
3249bd66eb0SPeter Brune   PetscFunctionBegin;
325f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
326f190f2fcSBarry Smith   if (func) linesearch->ops->precheck = func;
32786d74e61SPeter Brune   if (ctx) linesearch->precheckctx = ctx;
32886d74e61SPeter Brune   PetscFunctionReturn(0);
32986d74e61SPeter Brune }
33086d74e61SPeter Brune 
33186d74e61SPeter Brune /*@C
33253deb899SBarry Smith    SNESLineSearchGetPreCheck - Gets the pre-check function for the line search routine.
33386d74e61SPeter Brune 
334f899ff85SJose E. Roman    Input Parameter:
335f1c6b773SPeter Brune .  linesearch - the SNESLineSearch context
33686d74e61SPeter Brune 
33786d74e61SPeter Brune    Output Parameters:
33884238204SBarry Smith +  func       - [optional] function evaluation routine, see SNESLineSearchPreCheck() for calling sequence
339f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
34086d74e61SPeter Brune 
34186d74e61SPeter Brune    Level: intermediate
34286d74e61SPeter Brune 
34384238204SBarry Smith .seealso: SNESGetLineSearch(), SNESGetLineSearch(), SNESLineSearchPreCheck(), SNESLineSearchGetPostCheck(), SNESLineSearchSetPreCheck(), SNESLineSearchSetPostCheck()
34486d74e61SPeter Brune @*/
345f190f2fcSBarry Smith PetscErrorCode  SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void **ctx)
34686d74e61SPeter Brune {
3479bd66eb0SPeter Brune   PetscFunctionBegin;
348f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
349f190f2fcSBarry Smith   if (func) *func = linesearch->ops->precheck;
35086d74e61SPeter Brune   if (ctx) *ctx = linesearch->precheckctx;
35186d74e61SPeter Brune   PetscFunctionReturn(0);
35286d74e61SPeter Brune }
35386d74e61SPeter Brune 
35486d74e61SPeter Brune /*@C
355f190f2fcSBarry Smith    SNESLineSearchSetPostCheck - Sets a user function that is called after the line search has been applied to determine the step
356f190f2fcSBarry Smith        direction and length. Allows the user a chance to change or override the decision of the line search routine
35786d74e61SPeter Brune 
358f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
35986d74e61SPeter Brune 
36086d74e61SPeter Brune    Input Parameters:
361f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
36284238204SBarry Smith .  func - [optional] function evaluation routine, see SNESLineSearchPostCheck()  for the calling sequence
363f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
36486d74e61SPeter Brune 
36586d74e61SPeter Brune    Level: intermediate
36686d74e61SPeter Brune 
36784238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchPostCheck(), SNESLineSearchSetPreCheck(), SNESLineSearchGetPreCheck(), SNESLineSearchGetPostCheck()
36886d74e61SPeter Brune @*/
369f190f2fcSBarry Smith PetscErrorCode  SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void *ctx)
37086d74e61SPeter Brune {
37186d74e61SPeter Brune   PetscFunctionBegin;
372f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
373f190f2fcSBarry Smith   if (func) linesearch->ops->postcheck = func;
37486d74e61SPeter Brune   if (ctx) linesearch->postcheckctx = ctx;
37586d74e61SPeter Brune   PetscFunctionReturn(0);
37686d74e61SPeter Brune }
37786d74e61SPeter Brune 
37886d74e61SPeter Brune /*@C
379f1c6b773SPeter Brune    SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine.
38086d74e61SPeter Brune 
381f899ff85SJose E. Roman    Input Parameter:
382f1c6b773SPeter Brune .  linesearch - the SNESLineSearch context
38386d74e61SPeter Brune 
38486d74e61SPeter Brune    Output Parameters:
38584238204SBarry Smith +  func - [optional] function evaluation routine, see for the calling sequence SNESLineSearchPostCheck()
386f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
38786d74e61SPeter Brune 
38886d74e61SPeter Brune    Level: intermediate
38986d74e61SPeter Brune 
39084238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchPostCheck(), SNESLineSearchSetPreCheck()
39186d74e61SPeter Brune @*/
392f190f2fcSBarry Smith PetscErrorCode  SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void **ctx)
39386d74e61SPeter Brune {
3949bd66eb0SPeter Brune   PetscFunctionBegin;
395f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
396f190f2fcSBarry Smith   if (func) *func = linesearch->ops->postcheck;
39786d74e61SPeter Brune   if (ctx) *ctx = linesearch->postcheckctx;
39886d74e61SPeter Brune   PetscFunctionReturn(0);
39986d74e61SPeter Brune }
40086d74e61SPeter Brune 
401f40b411bSPeter Brune /*@
402f1c6b773SPeter Brune    SNESLineSearchPreCheck - Prepares the line search for being applied.
403f40b411bSPeter Brune 
404cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
405f40b411bSPeter Brune 
406f40b411bSPeter Brune    Input Parameters:
4077b1df9c1SPeter Brune +  linesearch - The linesearch instance.
4087b1df9c1SPeter Brune .  X - The current solution
4097b1df9c1SPeter Brune -  Y - The step direction
410f40b411bSPeter Brune 
411f40b411bSPeter Brune    Output Parameters:
4128e557f58SPeter Brune .  changed - Indicator that the precheck routine has changed anything
413f40b411bSPeter Brune 
414f190f2fcSBarry Smith    Level: developer
415f40b411bSPeter Brune 
41684238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchPostCheck(), SNESLineSearchSetPreCheck(), SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck()
417f40b411bSPeter Brune @*/
4187b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed)
419bf7f4e0aSPeter Brune {
420bf7f4e0aSPeter Brune   PetscFunctionBegin;
421bf7f4e0aSPeter Brune   *changed = PETSC_FALSE;
4226b2b7091SBarry Smith   if (linesearch->ops->precheck) {
4239566063dSJacob Faibussowitsch     PetscCall((*linesearch->ops->precheck)(linesearch, X, Y, changed, linesearch->precheckctx));
42438bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed,4);
425bf7f4e0aSPeter Brune   }
426bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
427bf7f4e0aSPeter Brune }
428bf7f4e0aSPeter Brune 
429f40b411bSPeter Brune /*@
430*ef46b1a6SStefano Zampini    SNESLineSearchPostCheck - Hook to modify step direction or updated solution after a successful linesearch
431f40b411bSPeter Brune 
432cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
433f40b411bSPeter Brune 
434f40b411bSPeter Brune    Input Parameters:
4357b1df9c1SPeter Brune +  linesearch - The linesearch context
4367b1df9c1SPeter Brune .  X - The last solution
4377b1df9c1SPeter Brune .  Y - The step direction
4387b1df9c1SPeter Brune -  W - The updated solution, W = X + lambda*Y for some lambda
439f40b411bSPeter Brune 
440f40b411bSPeter Brune    Output Parameters:
44178bcb3b5SPeter Brune +  changed_Y - Indicator if the direction Y has been changed.
44278bcb3b5SPeter Brune -  changed_W - Indicator if the new candidate solution W has been changed.
443f40b411bSPeter Brune 
444f190f2fcSBarry Smith    Level: developer
445f40b411bSPeter Brune 
44684238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck(), SNESLineSearchSetPrecheck(), SNESLineSearchGetPrecheck()
447f40b411bSPeter Brune @*/
4487b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch,Vec X,Vec Y,Vec W,PetscBool *changed_Y,PetscBool *changed_W)
449bf7f4e0aSPeter Brune {
450bf7f4e0aSPeter Brune   PetscFunctionBegin;
451bf7f4e0aSPeter Brune   *changed_Y = PETSC_FALSE;
452bf7f4e0aSPeter Brune   *changed_W = PETSC_FALSE;
4536b2b7091SBarry Smith   if (linesearch->ops->postcheck) {
4549566063dSJacob Faibussowitsch     PetscCall((*linesearch->ops->postcheck)(linesearch,X,Y,W,changed_Y,changed_W,linesearch->postcheckctx));
45538bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed_Y,5);
45638bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed_W,6);
45786d74e61SPeter Brune   }
45886d74e61SPeter Brune   PetscFunctionReturn(0);
45986d74e61SPeter Brune }
46086d74e61SPeter Brune 
46186d74e61SPeter Brune /*@C
46286d74e61SPeter Brune    SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration
46386d74e61SPeter Brune 
464cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
46586d74e61SPeter Brune 
4664165533cSJose E. Roman    Input Parameters:
46786d74e61SPeter Brune +  linesearch - linesearch context
46886d74e61SPeter Brune .  X - base state for this step
469907376e6SBarry Smith -  ctx - context for this function
47086d74e61SPeter Brune 
47197bb3fdcSJose E. Roman    Input/Output Parameter:
47297bb3fdcSJose E. Roman .  Y - correction, possibly modified
47397bb3fdcSJose E. Roman 
47497bb3fdcSJose E. Roman    Output Parameter:
47597bb3fdcSJose E. Roman .  changed - flag indicating that Y was modified
47686d74e61SPeter Brune 
47786d74e61SPeter Brune    Options Database Key:
478cd7522eaSPeter Brune +  -snes_linesearch_precheck_picard - activate this routine
479cd7522eaSPeter Brune -  -snes_linesearch_precheck_picard_angle - angle
48086d74e61SPeter Brune 
48186d74e61SPeter Brune    Level: advanced
48286d74e61SPeter Brune 
48386d74e61SPeter Brune    Notes:
48486d74e61SPeter Brune    This function should be passed to SNESLineSearchSetPreCheck()
48586d74e61SPeter Brune 
48686d74e61SPeter Brune    The justification for this method involves the linear convergence of a Picard iteration
48786d74e61SPeter Brune    so the Picard linearization should be provided in place of the "Jacobian". This correction
48886d74e61SPeter Brune    is generally not useful when using a Newton linearization.
48986d74e61SPeter Brune 
49086d74e61SPeter Brune    Reference:
49186d74e61SPeter Brune    Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology.
49286d74e61SPeter Brune 
49384238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchSetPreCheck()
49486d74e61SPeter Brune @*/
495f1c6b773SPeter Brune PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed,void *ctx)
49686d74e61SPeter Brune {
49786d74e61SPeter Brune   PetscReal      angle = *(PetscReal*)linesearch->precheckctx;
49886d74e61SPeter Brune   Vec            Ylast;
49986d74e61SPeter Brune   PetscScalar    dot;
50086d74e61SPeter Brune   PetscInt       iter;
50186d74e61SPeter Brune   PetscReal      ynorm,ylastnorm,theta,angle_radians;
50286d74e61SPeter Brune   SNES           snes;
50386d74e61SPeter Brune 
50486d74e61SPeter Brune   PetscFunctionBegin;
5059566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
5069566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject*)&Ylast));
50786d74e61SPeter Brune   if (!Ylast) {
5089566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(Y,&Ylast));
5099566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject)Ylast));
5109566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)Ylast));
51186d74e61SPeter Brune   }
5129566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes,&iter));
51386d74e61SPeter Brune   if (iter < 2) {
5149566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y,Ylast));
51586d74e61SPeter Brune     *changed = PETSC_FALSE;
51686d74e61SPeter Brune     PetscFunctionReturn(0);
51786d74e61SPeter Brune   }
51886d74e61SPeter Brune 
5199566063dSJacob Faibussowitsch   PetscCall(VecDot(Y,Ylast,&dot));
5209566063dSJacob Faibussowitsch   PetscCall(VecNorm(Y,NORM_2,&ynorm));
5219566063dSJacob Faibussowitsch   PetscCall(VecNorm(Ylast,NORM_2,&ylastnorm));
52286d74e61SPeter Brune   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
523255453a1SBarry Smith   theta         = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0));
52486d74e61SPeter Brune   angle_radians = angle * PETSC_PI / 180.;
52586d74e61SPeter Brune   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
52686d74e61SPeter Brune     /* Modify the step Y */
52786d74e61SPeter Brune     PetscReal alpha,ydiffnorm;
5289566063dSJacob Faibussowitsch     PetscCall(VecAXPY(Ylast,-1.0,Y));
5299566063dSJacob Faibussowitsch     PetscCall(VecNorm(Ylast,NORM_2,&ydiffnorm));
530f85e2ce2SBarry Smith     alpha = (ydiffnorm > .001*ylastnorm) ? ylastnorm / ydiffnorm : 1000.0;
5319566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y,Ylast));
5329566063dSJacob Faibussowitsch     PetscCall(VecScale(Y,alpha));
5339566063dSJacob 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));
534a47ec85fSBarry Smith     *changed = PETSC_TRUE;
53586d74e61SPeter Brune   } else {
5369566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes,"Angle %14.12e degrees exceeds threshold %14.12e, no correction applied\n",(double)(theta*180./PETSC_PI),(double)angle));
5379566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y,Ylast));
53886d74e61SPeter Brune     *changed = PETSC_FALSE;
539bf7f4e0aSPeter Brune   }
540bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
541bf7f4e0aSPeter Brune }
542bf7f4e0aSPeter Brune 
543f40b411bSPeter Brune /*@
544cd7522eaSPeter Brune    SNESLineSearchApply - Computes the line-search update.
545f40b411bSPeter Brune 
546f1c6b773SPeter Brune    Collective on SNESLineSearch
547f40b411bSPeter Brune 
548f40b411bSPeter Brune    Input Parameters:
5498e557f58SPeter Brune +  linesearch - The linesearch context
5508e557f58SPeter Brune -  Y - The search direction
551f40b411bSPeter Brune 
5526b867d5aSJose E. Roman    Input/Output Parameters:
5536b867d5aSJose E. Roman +  X - The current solution, on output the new solution
5546b867d5aSJose E. Roman .  F - The current function, on output the new function
5556b867d5aSJose E. Roman -  fnorm - The current norm, on output the new function norm
556f40b411bSPeter Brune 
557cd7522eaSPeter Brune    Options Database Keys:
558d4c6564cSPatrick Farrell + -snes_linesearch_type - basic, bt, l2, cp, nleqerr, shell
559dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches
5601fe24845SBarry Smith . -snes_linesearch_damping - The linesearch damping parameter, default is 1.0 (no damping)
5611fe24845SBarry Smith . -snes_linesearch_norms   - Turn on/off the linesearch norms computation (SNESLineSearchSetComputeNorms())
5623c7d6663SPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess
5633c7d6663SPeter Brune - -snes_linesearch_max_it - The number of iterations for iterative line searches
564cd7522eaSPeter Brune 
565cd7522eaSPeter Brune    Notes:
566cd7522eaSPeter Brune    This is typically called from within a SNESSolve() implementation in order to
567cd7522eaSPeter Brune    help with convergence of the nonlinear method.  Various SNES types use line searches
568cd7522eaSPeter Brune    in different ways, but the overarching theme is that a line search is used to determine
569cd7522eaSPeter Brune    an optimal damping parameter of a step at each iteration of the method.  Each
57084238204SBarry Smith    application of the line search may invoke SNESComputeFunction() several times, and
571cd7522eaSPeter Brune    therefore may be fairly expensive.
572cd7522eaSPeter Brune 
573f40b411bSPeter Brune    Level: Intermediate
574f40b411bSPeter Brune 
57584238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchCreate(), SNESLineSearchPreCheck(), SNESLineSearchPostCheck(), SNESSolve(), SNESComputeFunction(), SNESLineSearchSetComputeNorms(),
5761fe24845SBarry Smith           SNESLineSearchType, SNESLineSearchSetType()
577f40b411bSPeter Brune @*/
578bf388a1fSBarry Smith PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y)
579bf388a1fSBarry Smith {
580bf388a1fSBarry Smith   PetscFunctionBegin;
581f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
582bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
583bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
584064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(Y,VEC_CLASSID,5);
585bf7f4e0aSPeter Brune 
586422a814eSBarry Smith   linesearch->result = SNES_LINESEARCH_SUCCEEDED;
587bf7f4e0aSPeter Brune 
588bf7f4e0aSPeter Brune   linesearch->vec_sol    = X;
589bf7f4e0aSPeter Brune   linesearch->vec_update = Y;
590bf7f4e0aSPeter Brune   linesearch->vec_func   = F;
591bf7f4e0aSPeter Brune 
5929566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchSetUp(linesearch));
593bf7f4e0aSPeter Brune 
594f5af7f23SKarl Rupp   if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */
595bf7f4e0aSPeter Brune 
596f5af7f23SKarl Rupp   if (fnorm) linesearch->fnorm = *fnorm;
597f5af7f23SKarl Rupp   else {
5989566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &linesearch->fnorm));
599bf7f4e0aSPeter Brune   }
600bf7f4e0aSPeter Brune 
6019566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(SNESLINESEARCH_Apply,linesearch,X,F,Y));
602bf7f4e0aSPeter Brune 
6039566063dSJacob Faibussowitsch   PetscCall((*linesearch->ops->apply)(linesearch));
604bf7f4e0aSPeter Brune 
6059566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(SNESLINESEARCH_Apply,linesearch,X,F,Y));
606bf7f4e0aSPeter Brune 
607f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
608bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
609bf7f4e0aSPeter Brune }
610bf7f4e0aSPeter Brune 
611f40b411bSPeter Brune /*@
612f1c6b773SPeter Brune    SNESLineSearchDestroy - Destroys the line search instance.
613f40b411bSPeter Brune 
614f1c6b773SPeter Brune    Collective on SNESLineSearch
615f40b411bSPeter Brune 
616f40b411bSPeter Brune    Input Parameters:
6178e557f58SPeter Brune .  linesearch - The linesearch context
618f40b411bSPeter Brune 
61984238204SBarry Smith    Level: developer
620f40b411bSPeter Brune 
62184238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchCreate(), SNESLineSearchReset(), SNESDestroy()
622f40b411bSPeter Brune @*/
623bf388a1fSBarry Smith PetscErrorCode SNESLineSearchDestroy(SNESLineSearch * linesearch)
624bf388a1fSBarry Smith {
625bf7f4e0aSPeter Brune   PetscFunctionBegin;
626bf7f4e0aSPeter Brune   if (!*linesearch) PetscFunctionReturn(0);
627f1c6b773SPeter Brune   PetscValidHeaderSpecific((*linesearch),SNESLINESEARCH_CLASSID,1);
6289e5d0892SLisandro Dalcin   if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = NULL; PetscFunctionReturn(0);}
6299566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsViewOff((PetscObject)*linesearch));
6309566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchReset(*linesearch));
631f5af7f23SKarl Rupp   if ((*linesearch)->ops->destroy) (*linesearch)->ops->destroy(*linesearch);
6329566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&(*linesearch)->monitor));
6339566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchMonitorCancel((*linesearch)));
6349566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(linesearch));
635bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
636bf7f4e0aSPeter Brune }
637bf7f4e0aSPeter Brune 
638f40b411bSPeter Brune /*@
639dcf2fd19SBarry Smith    SNESLineSearchSetDefaultMonitor - Turns on/off printing useful information and debugging output about the line search.
640bf7f4e0aSPeter Brune 
641bf7f4e0aSPeter Brune    Input Parameters:
642dcf2fd19SBarry Smith +  linesearch - the linesearch object
643dcf2fd19SBarry Smith -  viewer - an ASCII PetscViewer or NULL to turn off monitor
644bf7f4e0aSPeter Brune 
645dcf2fd19SBarry Smith    Logically Collective on SNESLineSearch
646bf7f4e0aSPeter Brune 
647bf7f4e0aSPeter Brune    Options Database:
648dcf2fd19SBarry Smith .   -snes_linesearch_monitor [:filename] - enables the monitor
649bf7f4e0aSPeter Brune 
650bf7f4e0aSPeter Brune    Level: intermediate
651bf7f4e0aSPeter Brune 
652dcf2fd19SBarry Smith    Developer Note: This monitor is implemented differently than the other SNESLineSearchMonitors that are set with
653d12e167eSBarry Smith      SNESLineSearchMonitorSet() since it is called in many locations of the line search routines to display aspects of the
654d12e167eSBarry Smith      line search that are not visible to the other monitors.
655dcf2fd19SBarry Smith 
65684238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchGetDefaultMonitor(), PetscViewer, SNESLineSearchSetMonitor()
657bf7f4e0aSPeter Brune @*/
658dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchSetDefaultMonitor(SNESLineSearch linesearch, PetscViewer viewer)
659bf7f4e0aSPeter Brune {
660bf7f4e0aSPeter Brune   PetscFunctionBegin;
6619566063dSJacob Faibussowitsch   if (viewer) PetscCall(PetscObjectReference((PetscObject)viewer));
6629566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&linesearch->monitor));
663dcf2fd19SBarry Smith   linesearch->monitor = viewer;
664bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
665bf7f4e0aSPeter Brune }
666bf7f4e0aSPeter Brune 
667f40b411bSPeter Brune /*@
668dcf2fd19SBarry Smith    SNESLineSearchGetDefaultMonitor - Gets the PetscViewer instance for the line search monitor.
6696a388c36SPeter Brune 
670f190f2fcSBarry Smith    Input Parameter:
6718e557f58SPeter Brune .  linesearch - linesearch context
672f40b411bSPeter Brune 
673f190f2fcSBarry Smith    Output Parameter:
6748e557f58SPeter Brune .  monitor - monitor context
675f40b411bSPeter Brune 
676f40b411bSPeter Brune    Logically Collective on SNES
677f40b411bSPeter Brune 
678f40b411bSPeter Brune    Options Database Keys:
6798e557f58SPeter Brune .   -snes_linesearch_monitor - enables the monitor
680f40b411bSPeter Brune 
681f40b411bSPeter Brune    Level: intermediate
682f40b411bSPeter Brune 
68384238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchSetDefaultMonitor(), PetscViewer
684f40b411bSPeter Brune @*/
685dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchGetDefaultMonitor(SNESLineSearch linesearch, PetscViewer *monitor)
6866a388c36SPeter Brune {
6876a388c36SPeter Brune   PetscFunctionBegin;
688f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
6896a388c36SPeter Brune   *monitor = linesearch->monitor;
6906a388c36SPeter Brune   PetscFunctionReturn(0);
6916a388c36SPeter Brune }
6926a388c36SPeter Brune 
693dcf2fd19SBarry Smith /*@C
694dcf2fd19SBarry Smith    SNESLineSearchMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
695dcf2fd19SBarry Smith 
696dcf2fd19SBarry Smith    Collective on SNESLineSearch
697dcf2fd19SBarry Smith 
698dcf2fd19SBarry Smith    Input Parameters:
699dcf2fd19SBarry Smith +  ls - LineSearch object you wish to monitor
700dcf2fd19SBarry Smith .  name - the monitor type one is seeking
701dcf2fd19SBarry Smith .  help - message indicating what monitoring is done
702dcf2fd19SBarry Smith .  manual - manual page for the monitor
703dcf2fd19SBarry Smith .  monitor - the monitor function
704dcf2fd19SBarry 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
705dcf2fd19SBarry Smith 
706dcf2fd19SBarry Smith    Level: developer
707dcf2fd19SBarry Smith 
708dcf2fd19SBarry Smith .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
709dcf2fd19SBarry Smith           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
710dcf2fd19SBarry Smith           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
711d0609cedSBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHeadBegin(),
712dcf2fd19SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
713dcf2fd19SBarry Smith           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
714dcf2fd19SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
715dcf2fd19SBarry Smith @*/
716d12e167eSBarry Smith PetscErrorCode  SNESLineSearchMonitorSetFromOptions(SNESLineSearch ls,const char name[],const char help[], const char manual[],PetscErrorCode (*monitor)(SNESLineSearch,PetscViewerAndFormat*),PetscErrorCode (*monitorsetup)(SNESLineSearch,PetscViewerAndFormat*))
717dcf2fd19SBarry Smith {
718dcf2fd19SBarry Smith   PetscViewer       viewer;
719dcf2fd19SBarry Smith   PetscViewerFormat format;
720dcf2fd19SBarry Smith   PetscBool         flg;
721dcf2fd19SBarry Smith 
722dcf2fd19SBarry Smith   PetscFunctionBegin;
7239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)ls),((PetscObject) ls)->options,((PetscObject)ls)->prefix,name,&viewer,&format,&flg));
724dcf2fd19SBarry Smith   if (flg) {
725d12e167eSBarry Smith     PetscViewerAndFormat *vf;
7269566063dSJacob Faibussowitsch     PetscCall(PetscViewerAndFormatCreate(viewer,format,&vf));
7279566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)viewer));
728dcf2fd19SBarry Smith     if (monitorsetup) {
7299566063dSJacob Faibussowitsch       PetscCall((*monitorsetup)(ls,vf));
730dcf2fd19SBarry Smith     }
7319566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchMonitorSet(ls,(PetscErrorCode (*)(SNESLineSearch,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy));
732dcf2fd19SBarry Smith   }
733dcf2fd19SBarry Smith   PetscFunctionReturn(0);
734dcf2fd19SBarry Smith }
735dcf2fd19SBarry Smith 
736f40b411bSPeter Brune /*@
737f1c6b773SPeter Brune    SNESLineSearchSetFromOptions - Sets options for the line search
738f40b411bSPeter Brune 
739f899ff85SJose E. Roman    Input Parameter:
7408e557f58SPeter Brune .  linesearch - linesearch context
741f40b411bSPeter Brune 
742f40b411bSPeter Brune    Options Database Keys:
743d4c6564cSPatrick Farrell + -snes_linesearch_type <type> - basic, bt, l2, cp, nleqerr, shell
7443c7d6663SPeter Brune . -snes_linesearch_order <order> - 1, 2, 3.  Most types only support certain orders (bt supports 2 or 3)
7451fe24845SBarry Smith . -snes_linesearch_norms   - Turn on/off the linesearch norms for the basic linesearch typem (SNESLineSearchSetComputeNorms())
74671eef1aeSPeter Brune . -snes_linesearch_minlambda - The minimum step length
7471a9b3a06SPeter Brune . -snes_linesearch_maxstep - The maximum step size
7481a9b3a06SPeter Brune . -snes_linesearch_rtol - Relative tolerance for iterative line searches
7491a9b3a06SPeter Brune . -snes_linesearch_atol - Absolute tolerance for iterative line searches
7501a9b3a06SPeter Brune . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches
7511a9b3a06SPeter Brune . -snes_linesearch_max_it - The number of iterations for iterative line searches
752dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches
753dcf2fd19SBarry Smith . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine
7548e557f58SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter
755cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess.
7561a9b3a06SPeter Brune . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method
757d8d34be6SBarry Smith - -snes_linesearch_precheck_picard_angle - Angle used in Picard precheck method
758f40b411bSPeter Brune 
759f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
760f40b411bSPeter Brune 
761f40b411bSPeter Brune    Level: intermediate
762f40b411bSPeter Brune 
7631fe24845SBarry Smith .seealso: SNESLineSearchCreate(), SNESLineSearchSetOrder(), SNESLineSearchSetType(), SNESLineSearchSetTolerances(), SNESLineSearchSetDamping(), SNESLineSearchPreCheckPicard(),
7641fe24845SBarry Smith           SNESLineSearchType, SNESLineSearchSetComputeNorms()
765f40b411bSPeter Brune @*/
766bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch)
767bf388a1fSBarry Smith {
7681a4f838cSPeter Brune   const char        *deft = SNESLINESEARCHBASIC;
769bf7f4e0aSPeter Brune   char              type[256];
770bf7f4e0aSPeter Brune   PetscBool         flg, set;
771dcf2fd19SBarry Smith   PetscViewer       viewer;
772bf388a1fSBarry Smith 
773bf7f4e0aSPeter Brune   PetscFunctionBegin;
7749566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchRegisterAll());
775bf7f4e0aSPeter Brune 
776d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)linesearch);
777f5af7f23SKarl Rupp   if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name;
7789566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-snes_linesearch_type","Linesearch type","SNESLineSearchSetType",SNESLineSearchList,deft,type,256,&flg));
779bf7f4e0aSPeter Brune   if (flg) {
7809566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetType(linesearch,type));
781bf7f4e0aSPeter Brune   } else if (!((PetscObject)linesearch)->type_name) {
7829566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetType(linesearch,deft));
783bf7f4e0aSPeter Brune   }
784bf7f4e0aSPeter Brune 
7859566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)linesearch),((PetscObject) linesearch)->options,((PetscObject)linesearch)->prefix,"-snes_linesearch_monitor",&viewer,NULL,&set));
786dcf2fd19SBarry Smith   if (set) {
7879566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetDefaultMonitor(linesearch,viewer));
7889566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
789dcf2fd19SBarry Smith   }
7909566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchMonitorSetFromOptions(linesearch,"-snes_linesearch_monitor_solution_update","View correction at each iteration","SNESLineSearchMonitorSolutionUpdate",SNESLineSearchMonitorSolutionUpdate,NULL));
791bf7f4e0aSPeter Brune 
7921a9b3a06SPeter Brune   /* tolerances */
7939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_minlambda","Minimum step length","SNESLineSearchSetTolerances",linesearch->steptol,&linesearch->steptol,NULL));
7949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_maxstep","Maximum step size","SNESLineSearchSetTolerances",linesearch->maxstep,&linesearch->maxstep,NULL));
7959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_rtol","Relative tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->rtol,&linesearch->rtol,NULL));
7969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_atol","Absolute tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->atol,&linesearch->atol,NULL));
7979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_ltol","Change in lambda tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->ltol,&linesearch->ltol,NULL));
7989566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_linesearch_max_it","Maximum iterations for iterative line searches","SNESLineSearchSetTolerances",linesearch->max_its,&linesearch->max_its,NULL));
7997a35526eSPeter Brune 
8001a9b3a06SPeter Brune   /* damping parameters */
8019566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_damping","Line search damping and initial step guess","SNESLineSearchSetDamping",linesearch->damping,&linesearch->damping,NULL));
8021a9b3a06SPeter Brune 
8039566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_linesearch_keeplambda","Use previous lambda as damping","SNESLineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,NULL));
8041a9b3a06SPeter Brune 
8051a9b3a06SPeter Brune   /* precheck */
8069566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_linesearch_precheck_picard","Use a correction that sometimes improves convergence of Picard iteration","SNESLineSearchPreCheckPicard",flg,&flg,&set));
8077a35526eSPeter Brune   if (set) {
8087a35526eSPeter Brune     if (flg) {
8097a35526eSPeter Brune       linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */
810f5af7f23SKarl Rupp 
811d0609cedSBarry 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));
8129566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&linesearch->precheck_picard_angle));
8137a35526eSPeter Brune     } else {
8149566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchSetPreCheck(linesearch,NULL,NULL));
8157a35526eSPeter Brune     }
8167a35526eSPeter Brune   }
8179566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_linesearch_order","Order of approximation used in the line search","SNESLineSearchSetOrder",linesearch->order,&linesearch->order,NULL));
8189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_linesearch_norms","Compute final norms in line search","SNESLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,NULL));
8197a35526eSPeter Brune 
8205a9b6599SPeter Brune   if (linesearch->ops->setfromoptions) {
8219566063dSJacob Faibussowitsch     PetscCall((*linesearch->ops->setfromoptions)(PetscOptionsObject,linesearch));
8225a9b6599SPeter Brune   }
8235a9b6599SPeter Brune 
8249566063dSJacob Faibussowitsch   PetscCall(PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)linesearch));
825d0609cedSBarry Smith   PetscOptionsEnd();
826bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
827bf7f4e0aSPeter Brune }
828bf7f4e0aSPeter Brune 
829f40b411bSPeter Brune /*@
830f190f2fcSBarry Smith    SNESLineSearchView - Prints useful information about the line search
831f40b411bSPeter Brune 
832f40b411bSPeter Brune    Input Parameters:
8338e557f58SPeter Brune .  linesearch - linesearch context
834f40b411bSPeter Brune 
835f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
836f40b411bSPeter Brune 
837f40b411bSPeter Brune    Level: intermediate
838f40b411bSPeter Brune 
839dcf2fd19SBarry Smith .seealso: SNESLineSearchCreate()
840f40b411bSPeter Brune @*/
841bf388a1fSBarry Smith PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer)
842bf388a1fSBarry Smith {
8437f1410a3SPeter Brune   PetscBool      iascii;
844bf388a1fSBarry Smith 
845bf7f4e0aSPeter Brune   PetscFunctionBegin;
8467f1410a3SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
8477f1410a3SPeter Brune   if (!viewer) {
8489566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch),&viewer));
8497f1410a3SPeter Brune   }
8507f1410a3SPeter Brune   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
8517f1410a3SPeter Brune   PetscCheckSameComm(linesearch,1,viewer,2);
852f40b411bSPeter Brune 
8539566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
8547f1410a3SPeter Brune   if (iascii) {
8559566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)linesearch,viewer));
8567f1410a3SPeter Brune     if (linesearch->ops->view) {
8579566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
8589566063dSJacob Faibussowitsch       PetscCall((*linesearch->ops->view)(linesearch,viewer));
8599566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
8607f1410a3SPeter Brune     }
8619566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep,(double)linesearch->steptol));
8629566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol,(double)linesearch->atol,(double)linesearch->ltol));
8639566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D\n", linesearch->max_its));
8646b2b7091SBarry Smith     if (linesearch->ops->precheck) {
8656b2b7091SBarry Smith       if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) {
8669566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  using precheck step to speed up Picard convergence\n", linesearch->max_its));
8677f1410a3SPeter Brune       } else {
8689566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  using user-defined precheck step\n", linesearch->max_its));
8697f1410a3SPeter Brune       }
8707f1410a3SPeter Brune     }
8716b2b7091SBarry Smith     if (linesearch->ops->postcheck) {
8729566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  using user-defined postcheck step\n", linesearch->max_its));
8737f1410a3SPeter Brune     }
8747f1410a3SPeter Brune   }
875bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
876bf7f4e0aSPeter Brune }
877bf7f4e0aSPeter Brune 
878ea5d4fccSPeter Brune /*@C
879a80ff896SJed Brown    SNESLineSearchGetType - Gets the linesearch type
880a80ff896SJed Brown 
881a80ff896SJed Brown    Logically Collective on SNESLineSearch
882a80ff896SJed Brown 
883a80ff896SJed Brown    Input Parameters:
884a80ff896SJed Brown .  linesearch - linesearch context
885a80ff896SJed Brown 
886a80ff896SJed Brown    Output Parameters:
887a80ff896SJed Brown -  type - The type of line search, or NULL if not set
888a80ff896SJed Brown 
889a80ff896SJed Brown    Level: intermediate
890a80ff896SJed Brown 
891a80ff896SJed Brown .seealso: SNESLineSearchCreate(), SNESLineSearchType, SNESLineSearchSetFromOptions(), SNESLineSearchSetType()
892a80ff896SJed Brown @*/
893a80ff896SJed Brown PetscErrorCode SNESLineSearchGetType(SNESLineSearch linesearch, SNESLineSearchType *type)
894a80ff896SJed Brown {
895a80ff896SJed Brown   PetscFunctionBegin;
896a80ff896SJed Brown   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
897dadcf809SJacob Faibussowitsch   PetscValidPointer(type,2);
898a80ff896SJed Brown   *type = ((PetscObject)linesearch)->type_name;
899a80ff896SJed Brown   PetscFunctionReturn(0);
900a80ff896SJed Brown }
901a80ff896SJed Brown 
902a80ff896SJed Brown /*@C
903f1c6b773SPeter Brune    SNESLineSearchSetType - Sets the linesearch type
904f40b411bSPeter Brune 
905f190f2fcSBarry Smith    Logically Collective on SNESLineSearch
906f190f2fcSBarry Smith 
907f40b411bSPeter Brune    Input Parameters:
9088e557f58SPeter Brune +  linesearch - linesearch context
909f40b411bSPeter Brune -  type - The type of line search to be used
910f40b411bSPeter Brune 
911cd7522eaSPeter Brune    Available Types:
9121fe24845SBarry Smith +  SNESLINESEARCHBASIC - Simple damping line search, defaults to using the full Newton step
9131fe24845SBarry Smith .  SNESLINESEARCHBT - Backtracking line search over the L2 norm of the function
9141fe24845SBarry Smith .  SNESLINESEARCHL2 - Secant line search over the L2 norm of the function
9151fe24845SBarry Smith .  SNESLINESEARCHCP - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x)
9161fe24845SBarry Smith .  SNESLINESEARCHNLEQERR - Affine-covariant error-oriented linesearch
9171fe24845SBarry Smith -  SNESLINESEARCHSHELL - User provided SNESLineSearch implementation
9181fe24845SBarry Smith 
9191fe24845SBarry Smith    Options Database:
9201fe24845SBarry Smith .  -snes_linesearch_type <type> - basic, bt, l2, cp, nleqerr, shell
921cd7522eaSPeter Brune 
922f40b411bSPeter Brune    Level: intermediate
923f40b411bSPeter Brune 
924a80ff896SJed Brown .seealso: SNESLineSearchCreate(), SNESLineSearchType, SNESLineSearchSetFromOptions(), SNESLineSearchGetType()
925f40b411bSPeter Brune @*/
92619fd82e9SBarry Smith PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type)
927bf7f4e0aSPeter Brune {
928bf7f4e0aSPeter Brune   PetscBool      match;
9295f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(SNESLineSearch);
930bf7f4e0aSPeter Brune 
931bf7f4e0aSPeter Brune   PetscFunctionBegin;
932f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
933bf7f4e0aSPeter Brune   PetscValidCharPointer(type,2);
934bf7f4e0aSPeter Brune 
9359566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)linesearch,type,&match));
936bf7f4e0aSPeter Brune   if (match) PetscFunctionReturn(0);
937bf7f4e0aSPeter Brune 
9389566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(SNESLineSearchList,type,&r));
9395f80ce2aSJacob Faibussowitsch   PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type);
940bf7f4e0aSPeter Brune   /* Destroy the previous private linesearch context */
941bf7f4e0aSPeter Brune   if (linesearch->ops->destroy) {
9429566063dSJacob Faibussowitsch     PetscCall((*(linesearch)->ops->destroy)(linesearch));
9430298fd71SBarry Smith     linesearch->ops->destroy = NULL;
944bf7f4e0aSPeter Brune   }
945f1c6b773SPeter Brune   /* Reinitialize function pointers in SNESLineSearchOps structure */
9469e5d0892SLisandro Dalcin   linesearch->ops->apply          = NULL;
9479e5d0892SLisandro Dalcin   linesearch->ops->view           = NULL;
9489e5d0892SLisandro Dalcin   linesearch->ops->setfromoptions = NULL;
9499e5d0892SLisandro Dalcin   linesearch->ops->destroy        = NULL;
950bf7f4e0aSPeter Brune 
9519566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)linesearch,type));
9529566063dSJacob Faibussowitsch   PetscCall((*r)(linesearch));
953bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
954bf7f4e0aSPeter Brune }
955bf7f4e0aSPeter Brune 
956f40b411bSPeter Brune /*@
95778bcb3b5SPeter Brune    SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation.
958f40b411bSPeter Brune 
959f40b411bSPeter Brune    Input Parameters:
9608e557f58SPeter Brune +  linesearch - linesearch context
961f40b411bSPeter Brune -  snes - The snes instance
962f40b411bSPeter Brune 
96378bcb3b5SPeter Brune    Level: developer
96478bcb3b5SPeter Brune 
96578bcb3b5SPeter Brune    Notes:
966f190f2fcSBarry Smith    This happens automatically when the line search is obtained/created with
9677601faf0SJed Brown    SNESGetLineSearch().  This routine is therefore mainly called within SNES
96878bcb3b5SPeter Brune    implementations.
969f40b411bSPeter Brune 
9708141a3b9SPeter Brune    Level: developer
971f40b411bSPeter Brune 
972cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
973f40b411bSPeter Brune @*/
974bf388a1fSBarry Smith PetscErrorCode  SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes)
975bf388a1fSBarry Smith {
976bf7f4e0aSPeter Brune   PetscFunctionBegin;
977f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
978bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,2);
979bf7f4e0aSPeter Brune   linesearch->snes = snes;
980bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
981bf7f4e0aSPeter Brune }
982bf7f4e0aSPeter Brune 
983f40b411bSPeter Brune /*@
9848141a3b9SPeter Brune    SNESLineSearchGetSNES - Gets the SNES instance associated with the line search.
9858141a3b9SPeter Brune    Having an associated SNES is necessary because most line search implementations must be able to
9868141a3b9SPeter Brune    evaluate the function using SNESComputeFunction() for the associated SNES.  This routine
9878141a3b9SPeter Brune    is used in the line search implementations when one must get this associated SNES instance.
988f40b411bSPeter Brune 
989f40b411bSPeter Brune    Input Parameters:
9908e557f58SPeter Brune .  linesearch - linesearch context
991f40b411bSPeter Brune 
992f40b411bSPeter Brune    Output Parameters:
993f40b411bSPeter Brune .  snes - The snes instance
994f40b411bSPeter Brune 
9958141a3b9SPeter Brune    Level: developer
996f40b411bSPeter Brune 
997cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
998f40b411bSPeter Brune @*/
999bf388a1fSBarry Smith PetscErrorCode  SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes)
1000bf388a1fSBarry Smith {
1001bf7f4e0aSPeter Brune   PetscFunctionBegin;
1002f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
10036a388c36SPeter Brune   PetscValidPointer(snes,2);
1004bf7f4e0aSPeter Brune   *snes = linesearch->snes;
1005bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1006bf7f4e0aSPeter Brune }
1007bf7f4e0aSPeter Brune 
1008f40b411bSPeter Brune /*@
1009f1c6b773SPeter Brune    SNESLineSearchGetLambda - Gets the last linesearch steplength discovered.
1010f40b411bSPeter Brune 
1011f40b411bSPeter Brune    Input Parameters:
10128e557f58SPeter Brune .  linesearch - linesearch context
1013f40b411bSPeter Brune 
1014f40b411bSPeter Brune    Output Parameters:
1015cd7522eaSPeter Brune .  lambda - The last steplength computed during SNESLineSearchApply()
1016f40b411bSPeter Brune 
101778bcb3b5SPeter Brune    Level: advanced
101878bcb3b5SPeter Brune 
10198e557f58SPeter Brune    Notes:
10208e557f58SPeter Brune    This is useful in methods where the solver is ill-scaled and
102178bcb3b5SPeter Brune    requires some adaptive notion of the difference in scale between the
102278bcb3b5SPeter Brune    solution and the function.  For instance, SNESQN may be scaled by the
102378bcb3b5SPeter Brune    line search lambda using the argument -snes_qn_scaling ls.
102478bcb3b5SPeter Brune 
1025cd7522eaSPeter Brune .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetDamping(), SNESLineSearchApply()
1026f40b411bSPeter Brune @*/
1027f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda)
10286a388c36SPeter Brune {
10296a388c36SPeter Brune   PetscFunctionBegin;
1030f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1031534a8f05SLisandro Dalcin   PetscValidRealPointer(lambda, 2);
10326a388c36SPeter Brune   *lambda = linesearch->lambda;
10336a388c36SPeter Brune   PetscFunctionReturn(0);
10346a388c36SPeter Brune }
10356a388c36SPeter Brune 
1036f40b411bSPeter Brune /*@
1037f1c6b773SPeter Brune    SNESLineSearchSetLambda - Sets the linesearch steplength.
1038f40b411bSPeter Brune 
1039f40b411bSPeter Brune    Input Parameters:
10408e557f58SPeter Brune +  linesearch - linesearch context
1041f40b411bSPeter Brune -  lambda - The last steplength.
1042f40b411bSPeter Brune 
1043cd7522eaSPeter Brune    Notes:
1044f190f2fcSBarry Smith    This routine is typically used within implementations of SNESLineSearchApply()
1045cd7522eaSPeter Brune    to set the final steplength.  This routine (and SNESLineSearchGetLambda()) were
1046cd7522eaSPeter Brune    added in order to facilitate Quasi-Newton methods that use the previous steplength
1047cd7522eaSPeter Brune    as an inner scaling parameter.
1048cd7522eaSPeter Brune 
104978bcb3b5SPeter Brune    Level: advanced
1050f40b411bSPeter Brune 
1051f1c6b773SPeter Brune .seealso: SNESLineSearchGetLambda()
1052f40b411bSPeter Brune @*/
1053f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
10546a388c36SPeter Brune {
10556a388c36SPeter Brune   PetscFunctionBegin;
1056f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
10576a388c36SPeter Brune   linesearch->lambda = lambda;
10586a388c36SPeter Brune   PetscFunctionReturn(0);
10596a388c36SPeter Brune }
10606a388c36SPeter Brune 
1061f40b411bSPeter Brune /*@
10623c7d6663SPeter Brune    SNESLineSearchGetTolerances - Gets the tolerances for the linesearch.  These include
106378bcb3b5SPeter Brune    tolerances for the relative and absolute change in the function norm, the change
106478bcb3b5SPeter Brune    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
106578bcb3b5SPeter Brune    and the maximum number of iterations the line search procedure may take.
1066f40b411bSPeter Brune 
1067f899ff85SJose E. Roman    Input Parameter:
10688e557f58SPeter Brune .  linesearch - linesearch context
1069f40b411bSPeter Brune 
1070f40b411bSPeter Brune    Output Parameters:
1071516fe3c3SPeter Brune +  steptol - The minimum steplength
10726cc8e53bSPeter Brune .  maxstep - The maximum steplength
1073516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
1074516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
1075516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
1076516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
1077f40b411bSPeter Brune 
107878bcb3b5SPeter Brune    Level: intermediate
107978bcb3b5SPeter Brune 
108078bcb3b5SPeter Brune    Notes:
108178bcb3b5SPeter Brune    Different line searches may implement these parameters slightly differently as
10823c7d6663SPeter Brune    the type requires.
1083516fe3c3SPeter Brune 
1084f1c6b773SPeter Brune .seealso: SNESLineSearchSetTolerances()
1085f40b411bSPeter Brune @*/
1086f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
10876a388c36SPeter Brune {
10886a388c36SPeter Brune   PetscFunctionBegin;
1089f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1090516fe3c3SPeter Brune   if (steptol) {
1091534a8f05SLisandro Dalcin     PetscValidRealPointer(steptol, 2);
10926a388c36SPeter Brune     *steptol = linesearch->steptol;
1093516fe3c3SPeter Brune   }
1094516fe3c3SPeter Brune   if (maxstep) {
1095534a8f05SLisandro Dalcin     PetscValidRealPointer(maxstep, 3);
1096516fe3c3SPeter Brune     *maxstep = linesearch->maxstep;
1097516fe3c3SPeter Brune   }
1098516fe3c3SPeter Brune   if (rtol) {
1099534a8f05SLisandro Dalcin     PetscValidRealPointer(rtol, 4);
1100516fe3c3SPeter Brune     *rtol = linesearch->rtol;
1101516fe3c3SPeter Brune   }
1102516fe3c3SPeter Brune   if (atol) {
1103534a8f05SLisandro Dalcin     PetscValidRealPointer(atol, 5);
1104516fe3c3SPeter Brune     *atol = linesearch->atol;
1105516fe3c3SPeter Brune   }
1106516fe3c3SPeter Brune   if (ltol) {
1107534a8f05SLisandro Dalcin     PetscValidRealPointer(ltol, 6);
1108516fe3c3SPeter Brune     *ltol = linesearch->ltol;
1109516fe3c3SPeter Brune   }
1110516fe3c3SPeter Brune   if (max_its) {
1111534a8f05SLisandro Dalcin     PetscValidIntPointer(max_its, 7);
1112516fe3c3SPeter Brune     *max_its = linesearch->max_its;
1113516fe3c3SPeter Brune   }
11146a388c36SPeter Brune   PetscFunctionReturn(0);
11156a388c36SPeter Brune }
11166a388c36SPeter Brune 
1117f40b411bSPeter Brune /*@
11183c7d6663SPeter Brune    SNESLineSearchSetTolerances -  Gets the tolerances for the linesearch.  These include
111978bcb3b5SPeter Brune    tolerances for the relative and absolute change in the function norm, the change
112078bcb3b5SPeter Brune    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
112178bcb3b5SPeter Brune    and the maximum number of iterations the line search procedure may take.
1122f40b411bSPeter Brune 
1123f40b411bSPeter Brune    Input Parameters:
11248e557f58SPeter Brune +  linesearch - linesearch context
1125516fe3c3SPeter Brune .  steptol - The minimum steplength
11266cc8e53bSPeter Brune .  maxstep - The maximum steplength
1127516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
1128516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
1129516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
1130516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
1131f40b411bSPeter Brune 
113278bcb3b5SPeter Brune    Notes:
11333c7d6663SPeter Brune    The user may choose to not set any of the tolerances using PETSC_DEFAULT in
113478bcb3b5SPeter Brune    place of an argument.
1135f40b411bSPeter Brune 
113678bcb3b5SPeter Brune    Level: intermediate
1137516fe3c3SPeter Brune 
1138f1c6b773SPeter Brune .seealso: SNESLineSearchGetTolerances()
1139f40b411bSPeter Brune @*/
1140f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its)
11416a388c36SPeter Brune {
11426a388c36SPeter Brune   PetscFunctionBegin;
1143f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1144d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,steptol,2);
1145d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,maxstep,3);
1146d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,rtol,4);
1147d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,atol,5);
1148d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,ltol,6);
1149d3952184SSatish Balay   PetscValidLogicalCollectiveInt(linesearch,max_its,7);
1150d3952184SSatish Balay 
1151d3952184SSatish Balay   if (steptol!= PETSC_DEFAULT) {
11525f80ce2aSJacob Faibussowitsch     PetscCheck(steptol >= 0.0,PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Minimum step length %14.12e must be non-negative",(double)steptol);
11536a388c36SPeter Brune     linesearch->steptol = steptol;
1154d3952184SSatish Balay   }
1155d3952184SSatish Balay 
1156d3952184SSatish Balay   if (maxstep!= PETSC_DEFAULT) {
11575f80ce2aSJacob Faibussowitsch     PetscCheck(maxstep >= 0.0,PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum step length %14.12e must be non-negative",(double)maxstep);
1158516fe3c3SPeter Brune     linesearch->maxstep = maxstep;
1159d3952184SSatish Balay   }
1160d3952184SSatish Balay 
1161d3952184SSatish Balay   if (rtol != PETSC_DEFAULT) {
11625f80ce2aSJacob Faibussowitsch     PetscCheck(rtol >= 0.0 && 1.0 < rtol,PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Relative tolerance %14.12e must be non-negative and less than 1.0",(double)rtol);
1163516fe3c3SPeter Brune     linesearch->rtol = rtol;
1164d3952184SSatish Balay   }
1165d3952184SSatish Balay 
1166d3952184SSatish Balay   if (atol != PETSC_DEFAULT) {
11675f80ce2aSJacob Faibussowitsch     PetscCheck(atol >= 0.0,PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %14.12e must be non-negative",(double)atol);
1168516fe3c3SPeter Brune     linesearch->atol = atol;
1169d3952184SSatish Balay   }
1170d3952184SSatish Balay 
1171d3952184SSatish Balay   if (ltol != PETSC_DEFAULT) {
11725f80ce2aSJacob Faibussowitsch     PetscCheck(ltol >= 0.0,PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Lambda tolerance %14.12e must be non-negative",(double)ltol);
1173516fe3c3SPeter Brune     linesearch->ltol = ltol;
1174d3952184SSatish Balay   }
1175d3952184SSatish Balay 
1176d3952184SSatish Balay   if (max_its != PETSC_DEFAULT) {
11775f80ce2aSJacob Faibussowitsch     PetscCheck(max_its >= 0,PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",max_its);
1178516fe3c3SPeter Brune     linesearch->max_its = max_its;
1179d3952184SSatish Balay   }
11806a388c36SPeter Brune   PetscFunctionReturn(0);
11816a388c36SPeter Brune }
11826a388c36SPeter Brune 
1183f40b411bSPeter Brune /*@
1184f1c6b773SPeter Brune    SNESLineSearchGetDamping - Gets the line search damping parameter.
1185f40b411bSPeter Brune 
1186f40b411bSPeter Brune    Input Parameters:
11878e557f58SPeter Brune .  linesearch - linesearch context
1188f40b411bSPeter Brune 
1189f40b411bSPeter Brune    Output Parameters:
11908e557f58SPeter Brune .  damping - The damping parameter
1191f40b411bSPeter Brune 
119278bcb3b5SPeter Brune    Level: advanced
1193f40b411bSPeter Brune 
119478bcb3b5SPeter Brune .seealso: SNESLineSearchGetStepTolerance(), SNESQN
1195f40b411bSPeter Brune @*/
1196f40b411bSPeter Brune 
1197f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping)
11986a388c36SPeter Brune {
11996a388c36SPeter Brune   PetscFunctionBegin;
1200f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1201534a8f05SLisandro Dalcin   PetscValidRealPointer(damping, 2);
12026a388c36SPeter Brune   *damping = linesearch->damping;
12036a388c36SPeter Brune   PetscFunctionReturn(0);
12046a388c36SPeter Brune }
12056a388c36SPeter Brune 
1206f40b411bSPeter Brune /*@
1207fd292e60Sprj-    SNESLineSearchSetDamping - Sets the line search damping parameter.
1208f40b411bSPeter Brune 
1209f40b411bSPeter Brune    Input Parameters:
121003fd4120SBarry Smith +  linesearch - linesearch context
121103fd4120SBarry Smith -  damping - The damping parameter
1212f40b411bSPeter Brune 
121303fd4120SBarry Smith    Options Database:
121403fd4120SBarry Smith .   -snes_linesearch_damping
1215f40b411bSPeter Brune    Level: intermediate
1216f40b411bSPeter Brune 
1217cd7522eaSPeter Brune    Notes:
1218cd7522eaSPeter Brune    The basic line search merely takes the update step scaled by the damping parameter.
1219cd7522eaSPeter Brune    The use of the damping parameter in the l2 and cp line searches is much more subtle;
122078bcb3b5SPeter Brune    it is used as a starting point in calculating the secant step. However, the eventual
1221cd7522eaSPeter Brune    step may be of greater length than the damping parameter.  In the bt line search it is
1222cd7522eaSPeter Brune    used as the maximum possible step length, as the bt line search only backtracks.
1223cd7522eaSPeter Brune 
1224f1c6b773SPeter Brune .seealso: SNESLineSearchGetDamping()
1225f40b411bSPeter Brune @*/
1226f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping)
12276a388c36SPeter Brune {
12286a388c36SPeter Brune   PetscFunctionBegin;
1229f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
12306a388c36SPeter Brune   linesearch->damping = damping;
12316a388c36SPeter Brune   PetscFunctionReturn(0);
12326a388c36SPeter Brune }
12336a388c36SPeter Brune 
123459405d5eSPeter Brune /*@
123559405d5eSPeter Brune    SNESLineSearchGetOrder - Gets the line search approximation order.
123659405d5eSPeter Brune 
123759405d5eSPeter Brune    Input Parameters:
123878bcb3b5SPeter Brune .  linesearch - linesearch context
123959405d5eSPeter Brune 
124059405d5eSPeter Brune    Output Parameters:
12418e557f58SPeter Brune .  order - The order
124259405d5eSPeter Brune 
124378bcb3b5SPeter Brune    Possible Values for order:
12443c7d6663SPeter Brune +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
12453c7d6663SPeter Brune .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
12463c7d6663SPeter Brune -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
124778bcb3b5SPeter Brune 
124859405d5eSPeter Brune    Level: intermediate
124959405d5eSPeter Brune 
125059405d5eSPeter Brune .seealso: SNESLineSearchSetOrder()
125159405d5eSPeter Brune @*/
125259405d5eSPeter Brune 
1253b000cd8dSPeter Brune PetscErrorCode  SNESLineSearchGetOrder(SNESLineSearch linesearch,PetscInt *order)
125459405d5eSPeter Brune {
125559405d5eSPeter Brune   PetscFunctionBegin;
125659405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1257534a8f05SLisandro Dalcin   PetscValidIntPointer(order, 2);
125859405d5eSPeter Brune   *order = linesearch->order;
125959405d5eSPeter Brune   PetscFunctionReturn(0);
126059405d5eSPeter Brune }
126159405d5eSPeter Brune 
126259405d5eSPeter Brune /*@
12631f8196a2SJed Brown    SNESLineSearchSetOrder - Sets the maximum order of the polynomial fit used in the line search
126459405d5eSPeter Brune 
126559405d5eSPeter Brune    Input Parameters:
1266a2b725a8SWilliam Gropp +  linesearch - linesearch context
1267a2b725a8SWilliam Gropp -  order - The damping parameter
126859405d5eSPeter Brune 
126959405d5eSPeter Brune    Level: intermediate
127059405d5eSPeter Brune 
127178bcb3b5SPeter Brune    Possible Values for order:
12723c7d6663SPeter Brune +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
12733c7d6663SPeter Brune .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
12743c7d6663SPeter Brune -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
127578bcb3b5SPeter Brune 
127659405d5eSPeter Brune    Notes:
127759405d5eSPeter Brune    Variable orders are supported by the following line searches:
127878bcb3b5SPeter Brune +  bt - cubic and quadratic
127978bcb3b5SPeter Brune -  cp - linear and quadratic
128059405d5eSPeter Brune 
12811f8196a2SJed Brown .seealso: SNESLineSearchGetOrder(), SNESLineSearchSetDamping()
128259405d5eSPeter Brune @*/
1283b000cd8dSPeter Brune PetscErrorCode  SNESLineSearchSetOrder(SNESLineSearch linesearch,PetscInt order)
128459405d5eSPeter Brune {
128559405d5eSPeter Brune   PetscFunctionBegin;
128659405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
128759405d5eSPeter Brune   linesearch->order = order;
128859405d5eSPeter Brune   PetscFunctionReturn(0);
128959405d5eSPeter Brune }
129059405d5eSPeter Brune 
1291f40b411bSPeter Brune /*@
1292f1c6b773SPeter Brune    SNESLineSearchGetNorms - Gets the norms for for X, Y, and F.
1293f40b411bSPeter Brune 
1294f899ff85SJose E. Roman    Input Parameter:
129578bcb3b5SPeter Brune .  linesearch - linesearch context
1296f40b411bSPeter Brune 
1297f40b411bSPeter Brune    Output Parameters:
1298f40b411bSPeter Brune +  xnorm - The norm of the current solution
1299f40b411bSPeter Brune .  fnorm - The norm of the current function
1300f40b411bSPeter Brune -  ynorm - The norm of the current update
1301f40b411bSPeter Brune 
1302cd7522eaSPeter Brune    Notes:
1303cd7522eaSPeter Brune    This function is mainly called from SNES implementations.
1304cd7522eaSPeter Brune 
130578bcb3b5SPeter Brune    Level: developer
1306f40b411bSPeter Brune 
1307f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms() SNESLineSearchGetVecs()
1308f40b411bSPeter Brune @*/
1309f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm)
1310bf7f4e0aSPeter Brune {
1311bf7f4e0aSPeter Brune   PetscFunctionBegin;
1312f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1313f5af7f23SKarl Rupp   if (xnorm) *xnorm = linesearch->xnorm;
1314f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
1315f5af7f23SKarl Rupp   if (ynorm) *ynorm = linesearch->ynorm;
1316bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1317bf7f4e0aSPeter Brune }
1318bf7f4e0aSPeter Brune 
1319f40b411bSPeter Brune /*@
1320f1c6b773SPeter Brune    SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F.
1321f40b411bSPeter Brune 
1322f40b411bSPeter Brune    Input Parameters:
132378bcb3b5SPeter Brune +  linesearch - linesearch context
1324f40b411bSPeter Brune .  xnorm - The norm of the current solution
1325f40b411bSPeter Brune .  fnorm - The norm of the current function
1326f40b411bSPeter Brune -  ynorm - The norm of the current update
1327f40b411bSPeter Brune 
132878bcb3b5SPeter Brune    Level: advanced
1329f40b411bSPeter Brune 
1330f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1331f40b411bSPeter Brune @*/
1332f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
13336a388c36SPeter Brune {
13346a388c36SPeter Brune   PetscFunctionBegin;
1335f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
13366a388c36SPeter Brune   linesearch->xnorm = xnorm;
13376a388c36SPeter Brune   linesearch->fnorm = fnorm;
13386a388c36SPeter Brune   linesearch->ynorm = ynorm;
13396a388c36SPeter Brune   PetscFunctionReturn(0);
13406a388c36SPeter Brune }
13416a388c36SPeter Brune 
1342f40b411bSPeter Brune /*@
1343f1c6b773SPeter Brune    SNESLineSearchComputeNorms - Computes the norms of X, F, and Y.
1344f40b411bSPeter Brune 
1345f40b411bSPeter Brune    Input Parameters:
134678bcb3b5SPeter Brune .  linesearch - linesearch context
1347f40b411bSPeter Brune 
1348f40b411bSPeter Brune    Options Database Keys:
13498e557f58SPeter Brune .   -snes_linesearch_norms - turn norm computation on or off
1350f40b411bSPeter Brune 
1351f40b411bSPeter Brune    Level: intermediate
1352f40b411bSPeter Brune 
135378bcb3b5SPeter Brune .seealso: SNESLineSearchGetNorms, SNESLineSearchSetNorms(), SNESLineSearchSetComputeNorms()
1354f40b411bSPeter Brune @*/
1355f1c6b773SPeter Brune PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
13566a388c36SPeter Brune {
13579bd66eb0SPeter Brune   SNES           snes;
1358bf388a1fSBarry Smith 
13596a388c36SPeter Brune   PetscFunctionBegin;
13606a388c36SPeter Brune   if (linesearch->norms) {
13619bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
13629566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
13639566063dSJacob Faibussowitsch       PetscCall(VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
13649566063dSJacob Faibussowitsch       PetscCall(VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm));
13659566063dSJacob Faibussowitsch       PetscCall((*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm));
13669bd66eb0SPeter Brune     } else {
13679566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(linesearch->vec_func,   NORM_2, &linesearch->fnorm));
13689566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(linesearch->vec_sol,    NORM_2, &linesearch->xnorm));
13699566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm));
13709566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(linesearch->vec_func,     NORM_2, &linesearch->fnorm));
13719566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(linesearch->vec_sol,      NORM_2, &linesearch->xnorm));
13729566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(linesearch->vec_update,   NORM_2, &linesearch->ynorm));
13736a388c36SPeter Brune     }
13749bd66eb0SPeter Brune   }
13756a388c36SPeter Brune   PetscFunctionReturn(0);
13766a388c36SPeter Brune }
13776a388c36SPeter Brune 
13786f263ca3SPeter Brune /*@
13796f263ca3SPeter Brune    SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search.
13806f263ca3SPeter Brune 
13816f263ca3SPeter Brune    Input Parameters:
138278bcb3b5SPeter Brune +  linesearch  - linesearch context
138378bcb3b5SPeter Brune -  flg  - indicates whether or not to compute norms
13846f263ca3SPeter Brune 
13856f263ca3SPeter Brune    Options Database Keys:
13861f8196a2SJed Brown .   -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic linesearch
13876f263ca3SPeter Brune 
13886f263ca3SPeter Brune    Notes:
13891f8196a2SJed Brown    This is most relevant to the SNESLINESEARCHBASIC line search type since most line searches have a stopping criteria involving the norm.
13906f263ca3SPeter Brune 
13916f263ca3SPeter Brune    Level: intermediate
13926f263ca3SPeter Brune 
13931a4f838cSPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetNorms(), SNESLineSearchComputeNorms(), SNESLINESEARCHBASIC
13946f263ca3SPeter Brune @*/
13956f263ca3SPeter Brune PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg)
13966f263ca3SPeter Brune {
13976f263ca3SPeter Brune   PetscFunctionBegin;
13986f263ca3SPeter Brune   linesearch->norms = flg;
13996f263ca3SPeter Brune   PetscFunctionReturn(0);
14006f263ca3SPeter Brune }
14016f263ca3SPeter Brune 
1402f40b411bSPeter Brune /*@
1403f1c6b773SPeter Brune    SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context
1404f40b411bSPeter Brune 
1405f899ff85SJose E. Roman    Input Parameter:
140678bcb3b5SPeter Brune .  linesearch - linesearch context
1407f40b411bSPeter Brune 
1408f40b411bSPeter Brune    Output Parameters:
14096232e825SPeter Brune +  X - Solution vector
14106232e825SPeter Brune .  F - Function vector
14116232e825SPeter Brune .  Y - Search direction vector
14126232e825SPeter Brune .  W - Solution work vector
14136232e825SPeter Brune -  G - Function work vector
14146232e825SPeter Brune 
14157bba9028SPeter Brune    Notes:
14167bba9028SPeter Brune    At the beginning of a line search application, X should contain a
14176232e825SPeter Brune    solution and the vector F the function computed at X.  At the end of the
14186232e825SPeter Brune    line search application, X should contain the new solution, and F the
14196232e825SPeter Brune    function evaluated at the new solution.
1420f40b411bSPeter Brune 
14212a7a6963SBarry Smith    These vectors are owned by the SNESLineSearch and should not be destroyed by the caller
14222a7a6963SBarry Smith 
142378bcb3b5SPeter Brune    Level: advanced
1424f40b411bSPeter Brune 
1425f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1426f40b411bSPeter Brune @*/
1427bf388a1fSBarry Smith PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G)
1428bf388a1fSBarry Smith {
14296a388c36SPeter Brune   PetscFunctionBegin;
1430f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
14316a388c36SPeter Brune   if (X) {
14326a388c36SPeter Brune     PetscValidPointer(X, 2);
14336a388c36SPeter Brune     *X = linesearch->vec_sol;
14346a388c36SPeter Brune   }
14356a388c36SPeter Brune   if (F) {
14366a388c36SPeter Brune     PetscValidPointer(F, 3);
14376a388c36SPeter Brune     *F = linesearch->vec_func;
14386a388c36SPeter Brune   }
14396a388c36SPeter Brune   if (Y) {
14406a388c36SPeter Brune     PetscValidPointer(Y, 4);
14416a388c36SPeter Brune     *Y = linesearch->vec_update;
14426a388c36SPeter Brune   }
14436a388c36SPeter Brune   if (W) {
14446a388c36SPeter Brune     PetscValidPointer(W, 5);
14456a388c36SPeter Brune     *W = linesearch->vec_sol_new;
14466a388c36SPeter Brune   }
14476a388c36SPeter Brune   if (G) {
14486a388c36SPeter Brune     PetscValidPointer(G, 6);
14496a388c36SPeter Brune     *G = linesearch->vec_func_new;
14506a388c36SPeter Brune   }
14516a388c36SPeter Brune   PetscFunctionReturn(0);
14526a388c36SPeter Brune }
14536a388c36SPeter Brune 
1454f40b411bSPeter Brune /*@
1455f1c6b773SPeter Brune    SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context
1456f40b411bSPeter Brune 
1457f40b411bSPeter Brune    Input Parameters:
145878bcb3b5SPeter Brune +  linesearch - linesearch context
14596232e825SPeter Brune .  X - Solution vector
14606232e825SPeter Brune .  F - Function vector
14616232e825SPeter Brune .  Y - Search direction vector
14626232e825SPeter Brune .  W - Solution work vector
14636232e825SPeter Brune -  G - Function work vector
1464f40b411bSPeter Brune 
146578bcb3b5SPeter Brune    Level: advanced
1466f40b411bSPeter Brune 
1467f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms(), SNESLineSearchGetVecs()
1468f40b411bSPeter Brune @*/
1469bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G)
1470bf388a1fSBarry Smith {
14716a388c36SPeter Brune   PetscFunctionBegin;
1472f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
14736a388c36SPeter Brune   if (X) {
14746a388c36SPeter Brune     PetscValidHeaderSpecific(X,VEC_CLASSID,2);
14756a388c36SPeter Brune     linesearch->vec_sol = X;
14766a388c36SPeter Brune   }
14776a388c36SPeter Brune   if (F) {
14786a388c36SPeter Brune     PetscValidHeaderSpecific(F,VEC_CLASSID,3);
14796a388c36SPeter Brune     linesearch->vec_func = F;
14806a388c36SPeter Brune   }
14816a388c36SPeter Brune   if (Y) {
14826a388c36SPeter Brune     PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
14836a388c36SPeter Brune     linesearch->vec_update = Y;
14846a388c36SPeter Brune   }
14856a388c36SPeter Brune   if (W) {
14866a388c36SPeter Brune     PetscValidHeaderSpecific(W,VEC_CLASSID,5);
14876a388c36SPeter Brune     linesearch->vec_sol_new = W;
14886a388c36SPeter Brune   }
14896a388c36SPeter Brune   if (G) {
14906a388c36SPeter Brune     PetscValidHeaderSpecific(G,VEC_CLASSID,6);
14916a388c36SPeter Brune     linesearch->vec_func_new = G;
14926a388c36SPeter Brune   }
14936a388c36SPeter Brune   PetscFunctionReturn(0);
14946a388c36SPeter Brune }
14956a388c36SPeter Brune 
1496e7058c64SPeter Brune /*@C
1497f1c6b773SPeter Brune    SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1498e7058c64SPeter Brune    SNES options in the database.
1499e7058c64SPeter Brune 
1500cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
1501e7058c64SPeter Brune 
1502e7058c64SPeter Brune    Input Parameters:
1503e7058c64SPeter Brune +  snes - the SNES context
1504e7058c64SPeter Brune -  prefix - the prefix to prepend to all option names
1505e7058c64SPeter Brune 
1506e7058c64SPeter Brune    Notes:
1507e7058c64SPeter Brune    A hyphen (-) must NOT be given at the beginning of the prefix name.
1508e7058c64SPeter Brune    The first character of all runtime options is AUTOMATICALLY the hyphen.
1509e7058c64SPeter Brune 
1510e7058c64SPeter Brune    Level: advanced
1511e7058c64SPeter Brune 
1512e7058c64SPeter Brune .seealso: SNESGetOptionsPrefix()
1513e7058c64SPeter Brune @*/
1514f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[])
1515e7058c64SPeter Brune {
1516e7058c64SPeter Brune   PetscFunctionBegin;
1517f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
15189566063dSJacob Faibussowitsch   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix));
1519e7058c64SPeter Brune   PetscFunctionReturn(0);
1520e7058c64SPeter Brune }
1521e7058c64SPeter Brune 
1522e7058c64SPeter Brune /*@C
1523f1c6b773SPeter Brune    SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all
1524f1c6b773SPeter Brune    SNESLineSearch options in the database.
1525e7058c64SPeter Brune 
1526e7058c64SPeter Brune    Not Collective
1527e7058c64SPeter Brune 
1528e7058c64SPeter Brune    Input Parameter:
1529cd7522eaSPeter Brune .  linesearch - the SNESLineSearch context
1530e7058c64SPeter Brune 
1531e7058c64SPeter Brune    Output Parameter:
1532e7058c64SPeter Brune .  prefix - pointer to the prefix string used
1533e7058c64SPeter Brune 
15348e557f58SPeter Brune    Notes:
15358e557f58SPeter Brune    On the fortran side, the user should pass in a string 'prefix' of
1536e7058c64SPeter Brune    sufficient length to hold the prefix.
1537e7058c64SPeter Brune 
1538e7058c64SPeter Brune    Level: advanced
1539e7058c64SPeter Brune 
1540e7058c64SPeter Brune .seealso: SNESAppendOptionsPrefix()
1541e7058c64SPeter Brune @*/
1542f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[])
1543e7058c64SPeter Brune {
1544e7058c64SPeter Brune   PetscFunctionBegin;
1545f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
15469566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix));
1547e7058c64SPeter Brune   PetscFunctionReturn(0);
1548e7058c64SPeter Brune }
1549bf7f4e0aSPeter Brune 
15508d359177SBarry Smith /*@C
1551fa0ddf94SBarry Smith    SNESLineSearchSetWorkVecs - Gets work vectors for the line search.
1552f40b411bSPeter Brune 
1553d8d19677SJose E. Roman    Input Parameters:
1554f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
1555f40b411bSPeter Brune -  nwork - the number of work vectors
1556f40b411bSPeter Brune 
1557f40b411bSPeter Brune    Level: developer
1558f40b411bSPeter Brune 
1559fa0ddf94SBarry Smith .seealso: SNESSetWorkVecs()
1560f40b411bSPeter Brune @*/
1561fa0ddf94SBarry Smith PetscErrorCode  SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork)
1562bf7f4e0aSPeter Brune {
1563bf7f4e0aSPeter Brune   PetscFunctionBegin;
1564bf7f4e0aSPeter Brune   if (linesearch->vec_sol) {
15659566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work));
15668d359177SBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
1567bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1568bf7f4e0aSPeter Brune }
1569bf7f4e0aSPeter Brune 
1570f40b411bSPeter Brune /*@
1571422a814eSBarry Smith    SNESLineSearchGetReason - Gets the success/failure status of the last line search application
1572f40b411bSPeter Brune 
1573f40b411bSPeter Brune    Input Parameters:
157478bcb3b5SPeter Brune .  linesearch - linesearch context
1575f40b411bSPeter Brune 
1576f40b411bSPeter Brune    Output Parameters:
1577422a814eSBarry Smith .  result - The success or failure status
1578f40b411bSPeter Brune 
1579cd7522eaSPeter Brune    Notes:
1580c98378a5SBarry Smith    This is typically called after SNESLineSearchApply() in order to determine if the line-search failed
1581cd7522eaSPeter Brune    (and set the SNES convergence accordingly).
1582cd7522eaSPeter Brune 
1583f40b411bSPeter Brune    Level: intermediate
1584f40b411bSPeter Brune 
1585422a814eSBarry Smith .seealso: SNESLineSearchSetReason(), SNESLineSearchReason
1586f40b411bSPeter Brune @*/
1587422a814eSBarry Smith PetscErrorCode  SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result)
1588bf7f4e0aSPeter Brune {
1589bf7f4e0aSPeter Brune   PetscFunctionBegin;
1590f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1591422a814eSBarry Smith   PetscValidPointer(result, 2);
1592422a814eSBarry Smith   *result = linesearch->result;
1593bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1594bf7f4e0aSPeter Brune }
1595bf7f4e0aSPeter Brune 
1596f40b411bSPeter Brune /*@
1597422a814eSBarry Smith    SNESLineSearchSetReason - Sets the success/failure status of the last line search application
1598f40b411bSPeter Brune 
1599f40b411bSPeter Brune    Input Parameters:
160078bcb3b5SPeter Brune +  linesearch - linesearch context
1601422a814eSBarry Smith -  result - The success or failure status
1602f40b411bSPeter Brune 
1603cd7522eaSPeter Brune    Notes:
1604cd7522eaSPeter Brune    This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set
1605cd7522eaSPeter Brune    the success or failure of the line search method.
1606cd7522eaSPeter Brune 
160778bcb3b5SPeter Brune    Level: developer
1608f40b411bSPeter Brune 
1609422a814eSBarry Smith .seealso: SNESLineSearchGetSResult()
1610f40b411bSPeter Brune @*/
1611422a814eSBarry Smith PetscErrorCode  SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result)
16126a388c36SPeter Brune {
16136a388c36SPeter Brune   PetscFunctionBegin;
16145fd66863SKarl Rupp   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1615422a814eSBarry Smith   linesearch->result = result;
16166a388c36SPeter Brune   PetscFunctionReturn(0);
16176a388c36SPeter Brune }
16186a388c36SPeter Brune 
16199bd66eb0SPeter Brune /*@C
1620f1c6b773SPeter Brune    SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.
16219bd66eb0SPeter Brune 
16229bd66eb0SPeter Brune    Input Parameters:
16239bd66eb0SPeter Brune +  snes - nonlinear context obtained from SNESCreate()
16249bd66eb0SPeter Brune .  projectfunc - function for projecting the function to the bounds
16259bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
16269bd66eb0SPeter Brune 
16279bd66eb0SPeter Brune    Logically Collective on SNES
16289bd66eb0SPeter Brune 
16299bd66eb0SPeter Brune    Calling sequence of projectfunc:
16309bd66eb0SPeter Brune .vb
16319bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X)
16329bd66eb0SPeter Brune .ve
16339bd66eb0SPeter Brune 
16349bd66eb0SPeter Brune     Input parameters for projectfunc:
16359bd66eb0SPeter Brune +   snes - nonlinear context
16369bd66eb0SPeter Brune -   X - current solution
16379bd66eb0SPeter Brune 
1638cd7522eaSPeter Brune     Output parameters for projectfunc:
16399bd66eb0SPeter Brune .   X - Projected solution
16409bd66eb0SPeter Brune 
16419bd66eb0SPeter Brune    Calling sequence of normfunc:
16429bd66eb0SPeter Brune .vb
16439bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm)
16449bd66eb0SPeter Brune .ve
16459bd66eb0SPeter Brune 
1646cd7522eaSPeter Brune     Input parameters for normfunc:
16479bd66eb0SPeter Brune +   snes - nonlinear context
16489bd66eb0SPeter Brune .   X - current solution
16499bd66eb0SPeter Brune -   F - current residual
16509bd66eb0SPeter Brune 
1651cd7522eaSPeter Brune     Output parameters for normfunc:
16529bd66eb0SPeter Brune .   fnorm - VI-specific norm of the function
16539bd66eb0SPeter Brune 
1654cd7522eaSPeter Brune     Notes:
1655cd7522eaSPeter Brune     The VI solvers require projection of the solution to the feasible set.  projectfunc should implement this.
1656cd7522eaSPeter Brune 
1657cd7522eaSPeter Brune     The VI solvers require special evaluation of the function norm such that the norm is only calculated
1658cd7522eaSPeter Brune     on the inactive set.  This should be implemented by normfunc.
16599bd66eb0SPeter Brune 
16609bd66eb0SPeter Brune     Level: developer
16619bd66eb0SPeter Brune 
1662f1c6b773SPeter Brune .seealso: SNESLineSearchGetVIFunctions(), SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck()
16639bd66eb0SPeter Brune @*/
166425acbd8eSLisandro Dalcin PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc)
16659bd66eb0SPeter Brune {
16669bd66eb0SPeter Brune   PetscFunctionBegin;
1667f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
16689bd66eb0SPeter Brune   if (projectfunc) linesearch->ops->viproject = projectfunc;
16699bd66eb0SPeter Brune   if (normfunc) linesearch->ops->vinorm = normfunc;
16709bd66eb0SPeter Brune   PetscFunctionReturn(0);
16719bd66eb0SPeter Brune }
16729bd66eb0SPeter Brune 
16739bd66eb0SPeter Brune /*@C
1674f1c6b773SPeter Brune    SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.
16759bd66eb0SPeter Brune 
1676f899ff85SJose E. Roman    Input Parameter:
1677907376e6SBarry Smith .  linesearch - the line search context, obtain with SNESGetLineSearch()
16789bd66eb0SPeter Brune 
16799bd66eb0SPeter Brune    Output Parameters:
16809bd66eb0SPeter Brune +  projectfunc - function for projecting the function to the bounds
16819bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
16829bd66eb0SPeter Brune 
16839bd66eb0SPeter Brune    Logically Collective on SNES
16849bd66eb0SPeter Brune 
16859bd66eb0SPeter Brune     Level: developer
16869bd66eb0SPeter Brune 
1687f1c6b773SPeter Brune .seealso: SNESLineSearchSetVIFunctions(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
16889bd66eb0SPeter Brune @*/
168925acbd8eSLisandro Dalcin PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc)
16909bd66eb0SPeter Brune {
16919bd66eb0SPeter Brune   PetscFunctionBegin;
16929bd66eb0SPeter Brune   if (projectfunc) *projectfunc = linesearch->ops->viproject;
16939bd66eb0SPeter Brune   if (normfunc) *normfunc = linesearch->ops->vinorm;
16949bd66eb0SPeter Brune   PetscFunctionReturn(0);
16959bd66eb0SPeter Brune }
16969bd66eb0SPeter Brune 
1697bf7f4e0aSPeter Brune /*@C
16981c84c290SBarry Smith   SNESLineSearchRegister - See SNESLineSearchRegister()
1699bf7f4e0aSPeter Brune 
1700bf7f4e0aSPeter Brune   Level: advanced
1701bf7f4e0aSPeter Brune @*/
1702bdf89e91SBarry Smith PetscErrorCode  SNESLineSearchRegister(const char sname[],PetscErrorCode (*function)(SNESLineSearch))
1703bf7f4e0aSPeter Brune {
1704bf7f4e0aSPeter Brune   PetscFunctionBegin;
17059566063dSJacob Faibussowitsch   PetscCall(SNESInitializePackage());
17069566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&SNESLineSearchList,sname,function));
1707bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1708bf7f4e0aSPeter Brune }
1709