xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision a2b725a8db0d6bf6cc2a1c6df7dd8029aadfff6e)
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 
30dcf2fd19SBarry Smith .keywords: SNESLineSearch, nonlinear, set, monitor
31dcf2fd19SBarry Smith 
3284238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchMonitorDefault(), SNESLineSearchMonitorSet()
33dcf2fd19SBarry Smith @*/
34dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchMonitorCancel(SNESLineSearch ls)
35dcf2fd19SBarry Smith {
36dcf2fd19SBarry Smith   PetscErrorCode ierr;
37dcf2fd19SBarry Smith   PetscInt       i;
38dcf2fd19SBarry Smith 
39dcf2fd19SBarry Smith   PetscFunctionBegin;
40dcf2fd19SBarry Smith   PetscValidHeaderSpecific(ls,SNESLINESEARCH_CLASSID,1);
41dcf2fd19SBarry Smith   for (i=0; i<ls->numbermonitors; i++) {
42dcf2fd19SBarry Smith     if (ls->monitordestroy[i]) {
43dcf2fd19SBarry Smith       ierr = (*ls->monitordestroy[i])(&ls->monitorcontext[i]);CHKERRQ(ierr);
44dcf2fd19SBarry Smith     }
45dcf2fd19SBarry Smith   }
46dcf2fd19SBarry Smith   ls->numbermonitors = 0;
47dcf2fd19SBarry Smith   PetscFunctionReturn(0);
48dcf2fd19SBarry Smith }
49dcf2fd19SBarry Smith 
50dcf2fd19SBarry Smith /*@
51dcf2fd19SBarry Smith    SNESLineSearchMonitor - runs the user provided monitor routines, if they exist
52dcf2fd19SBarry Smith 
53dcf2fd19SBarry Smith    Collective on SNES
54dcf2fd19SBarry Smith 
55dcf2fd19SBarry Smith    Input Parameters:
56dcf2fd19SBarry Smith .  ls - the linesearch object
57dcf2fd19SBarry Smith 
58dcf2fd19SBarry Smith    Notes:
59dcf2fd19SBarry Smith    This routine is called by the SNES implementations.
60dcf2fd19SBarry Smith    It does not typically need to be called by the user.
61dcf2fd19SBarry Smith 
62dcf2fd19SBarry Smith    Level: developer
63dcf2fd19SBarry Smith 
6484238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchMonitorSet()
65dcf2fd19SBarry Smith @*/
66dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchMonitor(SNESLineSearch ls)
67dcf2fd19SBarry Smith {
68dcf2fd19SBarry Smith   PetscErrorCode ierr;
69dcf2fd19SBarry Smith   PetscInt       i,n = ls->numbermonitors;
70dcf2fd19SBarry Smith 
71dcf2fd19SBarry Smith   PetscFunctionBegin;
72dcf2fd19SBarry Smith   for (i=0; i<n; i++) {
73dcf2fd19SBarry Smith     ierr = (*ls->monitorftns[i])(ls,ls->monitorcontext[i]);CHKERRQ(ierr);
74dcf2fd19SBarry Smith   }
75dcf2fd19SBarry Smith   PetscFunctionReturn(0);
76dcf2fd19SBarry Smith }
77dcf2fd19SBarry Smith 
78dcf2fd19SBarry Smith /*@C
79dcf2fd19SBarry Smith    SNESLineSearchMonitorSet - Sets an ADDITIONAL function that is to be used at every
80dcf2fd19SBarry Smith    iteration of the nonlinear solver to display the iteration's
81dcf2fd19SBarry Smith    progress.
82dcf2fd19SBarry Smith 
83dcf2fd19SBarry Smith    Logically Collective on SNESLineSearch
84dcf2fd19SBarry Smith 
85dcf2fd19SBarry Smith    Input Parameters:
86dcf2fd19SBarry Smith +  ls - the SNESLineSearch context
87dcf2fd19SBarry Smith .  f - the monitor function
88dcf2fd19SBarry Smith .  mctx - [optional] user-defined context for private data for the
89dcf2fd19SBarry Smith           monitor routine (use NULL if no context is desired)
90dcf2fd19SBarry Smith -  monitordestroy - [optional] routine that frees monitor context
91dcf2fd19SBarry Smith           (may be NULL)
92dcf2fd19SBarry Smith 
93dcf2fd19SBarry Smith    Notes:
94dcf2fd19SBarry Smith    Several different monitoring routines may be set by calling
95dcf2fd19SBarry Smith    SNESLineSearchMonitorSet() multiple times; all will be called in the
96dcf2fd19SBarry Smith    order in which they were set.
97dcf2fd19SBarry Smith 
9895452b02SPatrick Sanan    Fortran Notes:
9995452b02SPatrick Sanan     Only a single monitor function can be set for each SNESLineSearch object
100dcf2fd19SBarry Smith 
101dcf2fd19SBarry Smith    Level: intermediate
102dcf2fd19SBarry Smith 
103dcf2fd19SBarry Smith .keywords: SNESLineSearch, nonlinear, set, monitor
104dcf2fd19SBarry Smith 
10584238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchMonitorDefault(), SNESLineSearchMonitorCancel()
106dcf2fd19SBarry Smith @*/
107dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchMonitorSet(SNESLineSearch ls,PetscErrorCode (*f)(SNESLineSearch,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
108dcf2fd19SBarry Smith {
10978064530SBarry Smith   PetscErrorCode ierr;
11078064530SBarry Smith   PetscInt       i;
11178064530SBarry Smith   PetscBool      identical;
11278064530SBarry Smith 
113dcf2fd19SBarry Smith   PetscFunctionBegin;
114dcf2fd19SBarry Smith   PetscValidHeaderSpecific(ls,SNESLINESEARCH_CLASSID,1);
11578064530SBarry Smith   for (i=0; i<ls->numbermonitors;i++) {
11678064530SBarry Smith     ierr = PetscMonitorCompare((PetscErrorCode (*)(void))f,mctx,monitordestroy,(PetscErrorCode (*)(void))ls->monitorftns[i],ls->monitorcontext[i],ls->monitordestroy[i],&identical);CHKERRQ(ierr);
11778064530SBarry Smith     if (identical) PetscFunctionReturn(0);
11878064530SBarry Smith   }
119dcf2fd19SBarry Smith   if (ls->numbermonitors >= MAXSNESLSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
120dcf2fd19SBarry Smith   ls->monitorftns[ls->numbermonitors]          = f;
121dcf2fd19SBarry Smith   ls->monitordestroy[ls->numbermonitors]   = monitordestroy;
122dcf2fd19SBarry Smith   ls->monitorcontext[ls->numbermonitors++] = (void*)mctx;
123dcf2fd19SBarry Smith   PetscFunctionReturn(0);
124dcf2fd19SBarry Smith }
125dcf2fd19SBarry Smith 
126dcf2fd19SBarry Smith /*@C
127dcf2fd19SBarry Smith    SNESLineSearchMonitorSolutionUpdate - Monitors each update a new function value the linesearch tries
128dcf2fd19SBarry Smith 
129dcf2fd19SBarry Smith    Collective on SNESLineSearch
130dcf2fd19SBarry Smith 
131dcf2fd19SBarry Smith    Input Parameters:
132dcf2fd19SBarry Smith +  ls - the SNES linesearch object
133d12e167eSBarry Smith -  vf - the context for the monitor, in this case it is an ASCII PetscViewer and format
134dcf2fd19SBarry Smith 
135dcf2fd19SBarry Smith    Level: intermediate
136dcf2fd19SBarry Smith 
137dcf2fd19SBarry Smith .keywords: SNES, nonlinear, default, monitor, norm
138dcf2fd19SBarry Smith 
13984238204SBarry Smith .seealso: SNESGetLineSearch(), SNESMonitorSet(), SNESMonitorSolution()
140dcf2fd19SBarry Smith @*/
141d12e167eSBarry Smith PetscErrorCode  SNESLineSearchMonitorSolutionUpdate(SNESLineSearch ls,PetscViewerAndFormat *vf)
142dcf2fd19SBarry Smith {
143dcf2fd19SBarry Smith   PetscErrorCode ierr;
144d12e167eSBarry Smith   PetscViewer    viewer = vf->viewer;
145dcf2fd19SBarry Smith   Vec            Y,W,G;
146dcf2fd19SBarry Smith 
147dcf2fd19SBarry Smith   PetscFunctionBegin;
148dcf2fd19SBarry Smith   ierr = SNESLineSearchGetVecs(ls,NULL,NULL,&Y,&W,&G);CHKERRQ(ierr);
149d12e167eSBarry Smith   ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr);
150dcf2fd19SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"LineSearch attempted update to solution \n");CHKERRQ(ierr);
151dcf2fd19SBarry Smith   ierr = VecView(Y,viewer);CHKERRQ(ierr);
152dcf2fd19SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"LineSearch attempted new solution \n");CHKERRQ(ierr);
153dcf2fd19SBarry Smith   ierr = VecView(W,viewer);CHKERRQ(ierr);
154dcf2fd19SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"LineSearch attempted updated function value\n");CHKERRQ(ierr);
155dcf2fd19SBarry Smith   ierr = VecView(G,viewer);CHKERRQ(ierr);
156d12e167eSBarry Smith   ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
157dcf2fd19SBarry Smith   PetscFunctionReturn(0);
158dcf2fd19SBarry Smith }
159dcf2fd19SBarry Smith 
160f40b411bSPeter Brune /*@
161cd7522eaSPeter Brune    SNESLineSearchCreate - Creates the line search context.
162f40b411bSPeter Brune 
163cd7522eaSPeter Brune    Logically Collective on Comm
164f40b411bSPeter Brune 
165f40b411bSPeter Brune    Input Parameters:
166cd7522eaSPeter Brune .  comm - MPI communicator for the line search (typically from the associated SNES context).
167f40b411bSPeter Brune 
168f40b411bSPeter Brune    Output Parameters:
1698e557f58SPeter Brune .  outlinesearch - the new linesearch context
170f40b411bSPeter Brune 
171162e0bf5SPeter Brune    Level: developer
172162e0bf5SPeter Brune 
173162e0bf5SPeter Brune    Notes:
174162e0bf5SPeter Brune    The preferred calling sequence for users is to use SNESGetLineSearch() to acquire the SNESLineSearch instance
175162e0bf5SPeter Brune    already associated with the SNES.  This function is for developer use.
176f40b411bSPeter Brune 
177cd7522eaSPeter Brune .keywords: LineSearch, create, context
178f40b411bSPeter Brune 
179162e0bf5SPeter Brune .seealso: LineSearchDestroy(), SNESGetLineSearch()
180f40b411bSPeter Brune @*/
181f40b411bSPeter Brune 
182bf388a1fSBarry Smith PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch)
183bf388a1fSBarry Smith {
184bf7f4e0aSPeter Brune   PetscErrorCode ierr;
185f1c6b773SPeter Brune   SNESLineSearch linesearch;
186bf388a1fSBarry Smith 
187bf7f4e0aSPeter Brune   PetscFunctionBegin;
188ea5d4fccSPeter Brune   PetscValidPointer(outlinesearch,2);
189f34a81feSMatthew G. Knepley   ierr = SNESInitializePackage();CHKERRQ(ierr);
1900298fd71SBarry Smith   *outlinesearch = NULL;
191f5af7f23SKarl Rupp 
19273107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(linesearch,SNESLINESEARCH_CLASSID, "SNESLineSearch","Linesearch","SNESLineSearch",comm,SNESLineSearchDestroy,SNESLineSearchView);CHKERRQ(ierr);
193bf7f4e0aSPeter Brune 
1940298fd71SBarry Smith   linesearch->vec_sol_new  = NULL;
1950298fd71SBarry Smith   linesearch->vec_func_new = NULL;
1960298fd71SBarry Smith   linesearch->vec_sol      = NULL;
1970298fd71SBarry Smith   linesearch->vec_func     = NULL;
1980298fd71SBarry Smith   linesearch->vec_update   = NULL;
1999bd66eb0SPeter Brune 
200bf7f4e0aSPeter Brune   linesearch->lambda       = 1.0;
201bf7f4e0aSPeter Brune   linesearch->fnorm        = 1.0;
202bf7f4e0aSPeter Brune   linesearch->ynorm        = 1.0;
203bf7f4e0aSPeter Brune   linesearch->xnorm        = 1.0;
204422a814eSBarry Smith   linesearch->result       = SNES_LINESEARCH_SUCCEEDED;
205bf7f4e0aSPeter Brune   linesearch->norms        = PETSC_TRUE;
206bf7f4e0aSPeter Brune   linesearch->keeplambda   = PETSC_FALSE;
207bf7f4e0aSPeter Brune   linesearch->damping      = 1.0;
208bf7f4e0aSPeter Brune   linesearch->maxstep      = 1e8;
209bf7f4e0aSPeter Brune   linesearch->steptol      = 1e-12;
210516fe3c3SPeter Brune   linesearch->rtol         = 1e-8;
211516fe3c3SPeter Brune   linesearch->atol         = 1e-15;
212516fe3c3SPeter Brune   linesearch->ltol         = 1e-8;
2130298fd71SBarry Smith   linesearch->precheckctx  = NULL;
2140298fd71SBarry Smith   linesearch->postcheckctx = NULL;
21559405d5eSPeter Brune   linesearch->max_its      = 1;
216bf7f4e0aSPeter Brune   linesearch->setupcalled  = PETSC_FALSE;
217bf7f4e0aSPeter Brune   *outlinesearch           = linesearch;
218bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
219bf7f4e0aSPeter Brune }
220bf7f4e0aSPeter Brune 
221f40b411bSPeter Brune /*@
22278bcb3b5SPeter Brune    SNESLineSearchSetUp - Prepares the line search for being applied by allocating
22378bcb3b5SPeter Brune    any required vectors.
224f40b411bSPeter Brune 
225cd7522eaSPeter Brune    Collective on SNESLineSearch
226f40b411bSPeter Brune 
227f40b411bSPeter Brune    Input Parameters:
228f40b411bSPeter Brune .  linesearch - The LineSearch instance.
229f40b411bSPeter Brune 
230cd7522eaSPeter Brune    Notes:
231f190f2fcSBarry Smith    For most cases, this needn't be called by users or outside of SNESLineSearchApply().
232cd7522eaSPeter Brune    The only current case where this is called outside of this is for the VI
23378bcb3b5SPeter Brune    solvers, which modify the solution and work vectors before the first call
234cd7522eaSPeter Brune    of SNESLineSearchApply, requiring the SNESLineSearch work vectors to be
235cd7522eaSPeter Brune    allocated upfront.
236cd7522eaSPeter Brune 
23778bcb3b5SPeter Brune    Level: advanced
238f40b411bSPeter Brune 
239f1c6b773SPeter Brune .keywords: SNESLineSearch, SetUp
240f40b411bSPeter Brune 
24184238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchReset()
242f40b411bSPeter Brune @*/
243f40b411bSPeter Brune 
244bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch)
245bf388a1fSBarry Smith {
246bf7f4e0aSPeter Brune   PetscErrorCode ierr;
247bf388a1fSBarry Smith 
248bf7f4e0aSPeter Brune   PetscFunctionBegin;
249bf7f4e0aSPeter Brune   if (!((PetscObject)linesearch)->type_name) {
2501a4f838cSPeter Brune     ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr);
251bf7f4e0aSPeter Brune   }
252bf7f4e0aSPeter Brune   if (!linesearch->setupcalled) {
2539bd66eb0SPeter Brune     if (!linesearch->vec_sol_new) {
254bf7f4e0aSPeter Brune       ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new);CHKERRQ(ierr);
2559bd66eb0SPeter Brune     }
2569bd66eb0SPeter Brune     if (!linesearch->vec_func_new) {
2579bd66eb0SPeter Brune       ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new);CHKERRQ(ierr);
2589bd66eb0SPeter Brune     }
259bf7f4e0aSPeter Brune     if (linesearch->ops->setup) {
260bf7f4e0aSPeter Brune       ierr = (*linesearch->ops->setup)(linesearch);CHKERRQ(ierr);
261bf7f4e0aSPeter Brune     }
262ed07d7d7SPeter Brune     if (!linesearch->ops->snesfunc) {ierr = SNESLineSearchSetFunction(linesearch,SNESComputeFunction);CHKERRQ(ierr);}
263bf7f4e0aSPeter Brune     linesearch->lambda      = linesearch->damping;
264bf7f4e0aSPeter Brune     linesearch->setupcalled = PETSC_TRUE;
265bf7f4e0aSPeter Brune   }
266bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
267bf7f4e0aSPeter Brune }
268bf7f4e0aSPeter Brune 
269f40b411bSPeter Brune 
270f40b411bSPeter Brune /*@
271f190f2fcSBarry Smith    SNESLineSearchReset - Undoes the SNESLineSearchSetUp() and deletes any Vecs or Mats allocated by the line search.
272f40b411bSPeter Brune 
273f1c6b773SPeter Brune    Collective on SNESLineSearch
274f40b411bSPeter Brune 
275f40b411bSPeter Brune    Input Parameters:
276f40b411bSPeter Brune .  linesearch - The LineSearch instance.
277f40b411bSPeter Brune 
27895452b02SPatrick Sanan    Notes:
27995452b02SPatrick Sanan     Usually only called by SNESReset()
280f190f2fcSBarry Smith 
281f190f2fcSBarry Smith    Level: developer
282f40b411bSPeter Brune 
283cd7522eaSPeter Brune .keywords: SNESLineSearch, Reset
284f40b411bSPeter Brune 
28584238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchSetUp()
286f40b411bSPeter Brune @*/
287f40b411bSPeter Brune 
288bf388a1fSBarry Smith PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch)
289bf388a1fSBarry Smith {
290bf7f4e0aSPeter Brune   PetscErrorCode ierr;
291bf388a1fSBarry Smith 
292bf7f4e0aSPeter Brune   PetscFunctionBegin;
293f5af7f23SKarl Rupp   if (linesearch->ops->reset) (*linesearch->ops->reset)(linesearch);
294f5af7f23SKarl Rupp 
295bf7f4e0aSPeter Brune   ierr = VecDestroy(&linesearch->vec_sol_new);CHKERRQ(ierr);
296bf7f4e0aSPeter Brune   ierr = VecDestroy(&linesearch->vec_func_new);CHKERRQ(ierr);
297bf7f4e0aSPeter Brune 
298bf7f4e0aSPeter Brune   ierr = VecDestroyVecs(linesearch->nwork, &linesearch->work);CHKERRQ(ierr);
299f5af7f23SKarl Rupp 
300bf7f4e0aSPeter Brune   linesearch->nwork       = 0;
301bf7f4e0aSPeter Brune   linesearch->setupcalled = PETSC_FALSE;
302bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
303bf7f4e0aSPeter Brune }
304bf7f4e0aSPeter Brune 
305ed07d7d7SPeter Brune /*@C
306f190f2fcSBarry Smith    SNESLineSearchSetFunction - Sets the function evaluation used by the SNES line search
307ed07d7d7SPeter Brune 
308ed07d7d7SPeter Brune    Input Parameters:
309ed07d7d7SPeter Brune .  linesearch - the SNESLineSearch context
310f190f2fcSBarry Smith +  func       - function evaluation routine
311ed07d7d7SPeter Brune 
312ed07d7d7SPeter Brune    Level: developer
313ed07d7d7SPeter Brune 
31495452b02SPatrick Sanan    Notes:
31595452b02SPatrick Sanan     This is used internally by PETSc and not called by users
316f190f2fcSBarry Smith 
317ed07d7d7SPeter Brune .keywords: get, linesearch, pre-check
318ed07d7d7SPeter Brune 
31984238204SBarry Smith .seealso: SNESGetLineSearch(), SNESSetFunction()
320ed07d7d7SPeter Brune @*/
321ed07d7d7SPeter Brune PetscErrorCode  SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES,Vec,Vec))
322ed07d7d7SPeter Brune {
323ed07d7d7SPeter Brune   PetscFunctionBegin;
324ed07d7d7SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
325ed07d7d7SPeter Brune   linesearch->ops->snesfunc = func;
326ed07d7d7SPeter Brune   PetscFunctionReturn(0);
327ed07d7d7SPeter Brune }
328ed07d7d7SPeter Brune 
32986d74e61SPeter Brune /*@C
330f190f2fcSBarry Smith    SNESLineSearchSetPreCheck - Sets a user function that is called after the initial search direction has been computed but
331df3898eeSBarry Smith          before the line search routine has been applied. Allows the user to adjust the result of (usually a linear solve) that
332f190f2fcSBarry Smith          determined the search direction.
33386d74e61SPeter Brune 
334f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
33586d74e61SPeter Brune 
33686d74e61SPeter Brune    Input Parameters:
337f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
33884238204SBarry Smith .  func - [optional] function evaluation routine, see SNESLineSearchPreCheck() for the 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 
34386d74e61SPeter Brune .keywords: set, linesearch, pre-check
34486d74e61SPeter Brune 
34584238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
34686d74e61SPeter Brune @*/
347f190f2fcSBarry Smith PetscErrorCode  SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void *ctx)
34886d74e61SPeter Brune {
3499bd66eb0SPeter Brune   PetscFunctionBegin;
350f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
351f190f2fcSBarry Smith   if (func) linesearch->ops->precheck = func;
35286d74e61SPeter Brune   if (ctx) linesearch->precheckctx = ctx;
35386d74e61SPeter Brune   PetscFunctionReturn(0);
35486d74e61SPeter Brune }
35586d74e61SPeter Brune 
35686d74e61SPeter Brune /*@C
35753deb899SBarry Smith    SNESLineSearchGetPreCheck - Gets the pre-check function for the line search routine.
35886d74e61SPeter Brune 
35986d74e61SPeter Brune    Input Parameters:
360f1c6b773SPeter Brune .  linesearch - the SNESLineSearch context
36186d74e61SPeter Brune 
36286d74e61SPeter Brune    Output Parameters:
36384238204SBarry Smith +  func       - [optional] function evaluation routine, see SNESLineSearchPreCheck() for calling sequence
364f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
36586d74e61SPeter Brune 
36686d74e61SPeter Brune    Level: intermediate
36786d74e61SPeter Brune 
36886d74e61SPeter Brune .keywords: get, linesearch, pre-check
36986d74e61SPeter Brune 
37084238204SBarry Smith .seealso: SNESGetLineSearch(), SNESGetLineSearch(), SNESLineSearchPreCheck(), SNESLineSearchGetPostCheck(), SNESLineSearchSetPreCheck(), SNESLineSearchSetPostCheck()
37186d74e61SPeter Brune @*/
372f190f2fcSBarry Smith PetscErrorCode  SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void **ctx)
37386d74e61SPeter Brune {
3749bd66eb0SPeter Brune   PetscFunctionBegin;
375f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
376f190f2fcSBarry Smith   if (func) *func = linesearch->ops->precheck;
37786d74e61SPeter Brune   if (ctx) *ctx = linesearch->precheckctx;
37886d74e61SPeter Brune   PetscFunctionReturn(0);
37986d74e61SPeter Brune }
38086d74e61SPeter Brune 
38186d74e61SPeter Brune /*@C
382f190f2fcSBarry Smith    SNESLineSearchSetPostCheck - Sets a user function that is called after the line search has been applied to determine the step
383f190f2fcSBarry Smith        direction and length. Allows the user a chance to change or override the decision of the line search routine
38486d74e61SPeter Brune 
385f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
38686d74e61SPeter Brune 
38786d74e61SPeter Brune    Input Parameters:
388f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
38984238204SBarry Smith .  func - [optional] function evaluation routine, see SNESLineSearchPostCheck()  for the calling sequence
390f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
39186d74e61SPeter Brune 
39286d74e61SPeter Brune    Level: intermediate
39386d74e61SPeter Brune 
39486d74e61SPeter Brune .keywords: set, linesearch, post-check
39586d74e61SPeter Brune 
39684238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchPostCheck(), SNESLineSearchSetPreCheck(), SNESLineSearchGetPreCheck(), SNESLineSearchGetPostCheck()
39786d74e61SPeter Brune @*/
398f190f2fcSBarry Smith PetscErrorCode  SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void *ctx)
39986d74e61SPeter Brune {
40086d74e61SPeter Brune   PetscFunctionBegin;
401f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
402f190f2fcSBarry Smith   if (func) linesearch->ops->postcheck = func;
40386d74e61SPeter Brune   if (ctx) linesearch->postcheckctx = ctx;
40486d74e61SPeter Brune   PetscFunctionReturn(0);
40586d74e61SPeter Brune }
40686d74e61SPeter Brune 
40786d74e61SPeter Brune /*@C
408f1c6b773SPeter Brune    SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine.
40986d74e61SPeter Brune 
41086d74e61SPeter Brune    Input Parameters:
411f1c6b773SPeter Brune .  linesearch - the SNESLineSearch context
41286d74e61SPeter Brune 
41386d74e61SPeter Brune    Output Parameters:
41484238204SBarry Smith +  func - [optional] function evaluation routine, see for the calling sequence SNESLineSearchPostCheck()
415f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
41686d74e61SPeter Brune 
41786d74e61SPeter Brune    Level: intermediate
41886d74e61SPeter Brune 
41986d74e61SPeter Brune .keywords: get, linesearch, post-check
42086d74e61SPeter Brune 
42184238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchPostCheck(), SNESLineSearchSetPreCheck()
42286d74e61SPeter Brune @*/
423f190f2fcSBarry Smith PetscErrorCode  SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void **ctx)
42486d74e61SPeter Brune {
4259bd66eb0SPeter Brune   PetscFunctionBegin;
426f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
427f190f2fcSBarry Smith   if (func) *func = linesearch->ops->postcheck;
42886d74e61SPeter Brune   if (ctx) *ctx = linesearch->postcheckctx;
42986d74e61SPeter Brune   PetscFunctionReturn(0);
43086d74e61SPeter Brune }
43186d74e61SPeter Brune 
432f40b411bSPeter Brune /*@
433f1c6b773SPeter Brune    SNESLineSearchPreCheck - Prepares the line search for being applied.
434f40b411bSPeter Brune 
435cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
436f40b411bSPeter Brune 
437f40b411bSPeter Brune    Input Parameters:
4387b1df9c1SPeter Brune +  linesearch - The linesearch instance.
4397b1df9c1SPeter Brune .  X - The current solution
4407b1df9c1SPeter Brune -  Y - The step direction
441f40b411bSPeter Brune 
442f40b411bSPeter Brune    Output Parameters:
4438e557f58SPeter Brune .  changed - Indicator that the precheck routine has changed anything
444f40b411bSPeter Brune 
445f190f2fcSBarry Smith    Level: developer
446f40b411bSPeter Brune 
447f1c6b773SPeter Brune .keywords: SNESLineSearch, Create
448f40b411bSPeter Brune 
44984238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchPostCheck(), SNESLineSearchSetPreCheck(), SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck()
450f40b411bSPeter Brune @*/
4517b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed)
452bf7f4e0aSPeter Brune {
453bf7f4e0aSPeter Brune   PetscErrorCode ierr;
4545fd66863SKarl Rupp 
455bf7f4e0aSPeter Brune   PetscFunctionBegin;
456bf7f4e0aSPeter Brune   *changed = PETSC_FALSE;
4576b2b7091SBarry Smith   if (linesearch->ops->precheck) {
4586b2b7091SBarry Smith     ierr = (*linesearch->ops->precheck)(linesearch, X, Y, changed, linesearch->precheckctx);CHKERRQ(ierr);
45938bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed,4);
460bf7f4e0aSPeter Brune   }
461bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
462bf7f4e0aSPeter Brune }
463bf7f4e0aSPeter Brune 
464f40b411bSPeter Brune /*@
465f1c6b773SPeter Brune    SNESLineSearchPostCheck - Prepares the line search for being applied.
466f40b411bSPeter Brune 
467cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
468f40b411bSPeter Brune 
469f40b411bSPeter Brune    Input Parameters:
4707b1df9c1SPeter Brune +  linesearch - The linesearch context
4717b1df9c1SPeter Brune .  X - The last solution
4727b1df9c1SPeter Brune .  Y - The step direction
4737b1df9c1SPeter Brune -  W - The updated solution, W = X + lambda*Y for some lambda
474f40b411bSPeter Brune 
475f40b411bSPeter Brune    Output Parameters:
47678bcb3b5SPeter Brune +  changed_Y - Indicator if the direction Y has been changed.
47778bcb3b5SPeter Brune -  changed_W - Indicator if the new candidate solution W has been changed.
478f40b411bSPeter Brune 
479f190f2fcSBarry Smith    Level: developer
480f40b411bSPeter Brune 
481f1c6b773SPeter Brune .keywords: SNESLineSearch, Create
482f40b411bSPeter Brune 
48384238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck(), SNESLineSearchSetPrecheck(), SNESLineSearchGetPrecheck()
484f40b411bSPeter Brune @*/
4857b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch,Vec X,Vec Y,Vec W,PetscBool *changed_Y,PetscBool *changed_W)
486bf7f4e0aSPeter Brune {
487bf7f4e0aSPeter Brune   PetscErrorCode ierr;
488bf388a1fSBarry Smith 
489bf7f4e0aSPeter Brune   PetscFunctionBegin;
490bf7f4e0aSPeter Brune   *changed_Y = PETSC_FALSE;
491bf7f4e0aSPeter Brune   *changed_W = PETSC_FALSE;
4926b2b7091SBarry Smith   if (linesearch->ops->postcheck) {
4936b2b7091SBarry Smith     ierr = (*linesearch->ops->postcheck)(linesearch,X,Y,W,changed_Y,changed_W,linesearch->postcheckctx);CHKERRQ(ierr);
49438bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed_Y,5);
49538bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed_W,6);
49686d74e61SPeter Brune   }
49786d74e61SPeter Brune   PetscFunctionReturn(0);
49886d74e61SPeter Brune }
49986d74e61SPeter Brune 
50086d74e61SPeter Brune /*@C
50186d74e61SPeter Brune    SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration
50286d74e61SPeter Brune 
503cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
50486d74e61SPeter Brune 
50586d74e61SPeter Brune    Input Arguments:
50686d74e61SPeter Brune +  linesearch - linesearch context
50786d74e61SPeter Brune .  X - base state for this step
50886d74e61SPeter Brune .  Y - initial correction
509907376e6SBarry Smith -  ctx - context for this function
51086d74e61SPeter Brune 
51186d74e61SPeter Brune    Output Arguments:
51286d74e61SPeter Brune +  Y - correction, possibly modified
51386d74e61SPeter Brune -  changed - flag indicating that Y was modified
51486d74e61SPeter Brune 
51586d74e61SPeter Brune    Options Database Key:
516cd7522eaSPeter Brune +  -snes_linesearch_precheck_picard - activate this routine
517cd7522eaSPeter Brune -  -snes_linesearch_precheck_picard_angle - angle
51886d74e61SPeter Brune 
51986d74e61SPeter Brune    Level: advanced
52086d74e61SPeter Brune 
52186d74e61SPeter Brune    Notes:
52286d74e61SPeter Brune    This function should be passed to SNESLineSearchSetPreCheck()
52386d74e61SPeter Brune 
52486d74e61SPeter Brune    The justification for this method involves the linear convergence of a Picard iteration
52586d74e61SPeter Brune    so the Picard linearization should be provided in place of the "Jacobian". This correction
52686d74e61SPeter Brune    is generally not useful when using a Newton linearization.
52786d74e61SPeter Brune 
52886d74e61SPeter Brune    Reference:
52986d74e61SPeter Brune    Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology.
53086d74e61SPeter Brune 
53184238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchSetPreCheck()
53286d74e61SPeter Brune @*/
533f1c6b773SPeter Brune PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed,void *ctx)
53486d74e61SPeter Brune {
53586d74e61SPeter Brune   PetscErrorCode ierr;
53686d74e61SPeter Brune   PetscReal      angle = *(PetscReal*)linesearch->precheckctx;
53786d74e61SPeter Brune   Vec            Ylast;
53886d74e61SPeter Brune   PetscScalar    dot;
53986d74e61SPeter Brune   PetscInt       iter;
54086d74e61SPeter Brune   PetscReal      ynorm,ylastnorm,theta,angle_radians;
54186d74e61SPeter Brune   SNES           snes;
54286d74e61SPeter Brune 
54386d74e61SPeter Brune   PetscFunctionBegin;
544f1c6b773SPeter Brune   ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
54586d74e61SPeter Brune   ierr = PetscObjectQuery((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject*)&Ylast);CHKERRQ(ierr);
54686d74e61SPeter Brune   if (!Ylast) {
54786d74e61SPeter Brune     ierr = VecDuplicate(Y,&Ylast);CHKERRQ(ierr);
54886d74e61SPeter Brune     ierr = PetscObjectCompose((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject)Ylast);CHKERRQ(ierr);
54986d74e61SPeter Brune     ierr = PetscObjectDereference((PetscObject)Ylast);CHKERRQ(ierr);
55086d74e61SPeter Brune   }
55186d74e61SPeter Brune   ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr);
55286d74e61SPeter Brune   if (iter < 2) {
55386d74e61SPeter Brune     ierr     = VecCopy(Y,Ylast);CHKERRQ(ierr);
55486d74e61SPeter Brune     *changed = PETSC_FALSE;
55586d74e61SPeter Brune     PetscFunctionReturn(0);
55686d74e61SPeter Brune   }
55786d74e61SPeter Brune 
55886d74e61SPeter Brune   ierr = VecDot(Y,Ylast,&dot);CHKERRQ(ierr);
55986d74e61SPeter Brune   ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr);
56086d74e61SPeter Brune   ierr = VecNorm(Ylast,NORM_2,&ylastnorm);CHKERRQ(ierr);
56186d74e61SPeter Brune   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
562255453a1SBarry Smith   theta         = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0));
56386d74e61SPeter Brune   angle_radians = angle * PETSC_PI / 180.;
56486d74e61SPeter Brune   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
56586d74e61SPeter Brune     /* Modify the step Y */
56686d74e61SPeter Brune     PetscReal alpha,ydiffnorm;
56786d74e61SPeter Brune     ierr  = VecAXPY(Ylast,-1.0,Y);CHKERRQ(ierr);
56886d74e61SPeter Brune     ierr  = VecNorm(Ylast,NORM_2,&ydiffnorm);CHKERRQ(ierr);
569f85e2ce2SBarry Smith     alpha = (ydiffnorm > .001*ylastnorm) ? ylastnorm / ydiffnorm : 1000.0;
57086d74e61SPeter Brune     ierr  = VecCopy(Y,Ylast);CHKERRQ(ierr);
57186d74e61SPeter Brune     ierr  = VecScale(Y,alpha);CHKERRQ(ierr);
572c69d1a72SBarry Smith     ierr  = PetscInfo3(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);CHKERRQ(ierr);
573a47ec85fSBarry Smith     *changed = PETSC_TRUE;
57486d74e61SPeter Brune   } else {
575c69d1a72SBarry Smith     ierr     = PetscInfo2(snes,"Angle %14.12e degrees exceeds threshold %14.12e, no correction applied\n",(double)(theta*180./PETSC_PI),(double)angle);CHKERRQ(ierr);
57686d74e61SPeter Brune     ierr     = VecCopy(Y,Ylast);CHKERRQ(ierr);
57786d74e61SPeter Brune     *changed = PETSC_FALSE;
578bf7f4e0aSPeter Brune   }
579bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
580bf7f4e0aSPeter Brune }
581bf7f4e0aSPeter Brune 
582f40b411bSPeter Brune /*@
583cd7522eaSPeter Brune    SNESLineSearchApply - Computes the line-search update.
584f40b411bSPeter Brune 
585f1c6b773SPeter Brune    Collective on SNESLineSearch
586f40b411bSPeter Brune 
587f40b411bSPeter Brune    Input Parameters:
5888e557f58SPeter Brune +  linesearch - The linesearch context
5898e557f58SPeter Brune .  X - The current solution
5908e557f58SPeter Brune .  F - The current function
5918e557f58SPeter Brune .  fnorm - The current norm
5928e557f58SPeter Brune -  Y - The search direction
593f40b411bSPeter Brune 
594f40b411bSPeter Brune    Output Parameters:
5958e557f58SPeter Brune +  X - The new solution
5968e557f58SPeter Brune .  F - The new function
5978e557f58SPeter Brune -  fnorm - The new function norm
598f40b411bSPeter Brune 
599cd7522eaSPeter Brune    Options Database Keys:
600d4c6564cSPatrick Farrell + -snes_linesearch_type - basic, bt, l2, cp, nleqerr, shell
601dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches
6021fe24845SBarry Smith . -snes_linesearch_damping - The linesearch damping parameter, default is 1.0 (no damping)
6031fe24845SBarry Smith . -snes_linesearch_norms   - Turn on/off the linesearch norms computation (SNESLineSearchSetComputeNorms())
6043c7d6663SPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess
6053c7d6663SPeter Brune - -snes_linesearch_max_it - The number of iterations for iterative line searches
606cd7522eaSPeter Brune 
607cd7522eaSPeter Brune    Notes:
608cd7522eaSPeter Brune    This is typically called from within a SNESSolve() implementation in order to
609cd7522eaSPeter Brune    help with convergence of the nonlinear method.  Various SNES types use line searches
610cd7522eaSPeter Brune    in different ways, but the overarching theme is that a line search is used to determine
611cd7522eaSPeter Brune    an optimal damping parameter of a step at each iteration of the method.  Each
61284238204SBarry Smith    application of the line search may invoke SNESComputeFunction() several times, and
613cd7522eaSPeter Brune    therefore may be fairly expensive.
614cd7522eaSPeter Brune 
615f40b411bSPeter Brune    Level: Intermediate
616f40b411bSPeter Brune 
617f1c6b773SPeter Brune .keywords: SNESLineSearch, Create
618f40b411bSPeter Brune 
61984238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchCreate(), SNESLineSearchPreCheck(), SNESLineSearchPostCheck(), SNESSolve(), SNESComputeFunction(), SNESLineSearchSetComputeNorms(),
6201fe24845SBarry Smith           SNESLineSearchType, SNESLineSearchSetType()
621f40b411bSPeter Brune @*/
622bf388a1fSBarry Smith PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y)
623bf388a1fSBarry Smith {
624bf7f4e0aSPeter Brune   PetscErrorCode ierr;
625bf7f4e0aSPeter Brune 
626bf388a1fSBarry Smith   PetscFunctionBegin;
627f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
628bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
629bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
630bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
631bf7f4e0aSPeter Brune 
632422a814eSBarry Smith   linesearch->result = SNES_LINESEARCH_SUCCEEDED;
633bf7f4e0aSPeter Brune 
634bf7f4e0aSPeter Brune   linesearch->vec_sol    = X;
635bf7f4e0aSPeter Brune   linesearch->vec_update = Y;
636bf7f4e0aSPeter Brune   linesearch->vec_func   = F;
637bf7f4e0aSPeter Brune 
638f1c6b773SPeter Brune   ierr = SNESLineSearchSetUp(linesearch);CHKERRQ(ierr);
639bf7f4e0aSPeter Brune 
640f5af7f23SKarl Rupp   if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */
641bf7f4e0aSPeter Brune 
642f5af7f23SKarl Rupp   if (fnorm) linesearch->fnorm = *fnorm;
643f5af7f23SKarl Rupp   else {
644bf7f4e0aSPeter Brune     ierr = VecNorm(F, NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
645bf7f4e0aSPeter Brune   }
646bf7f4e0aSPeter Brune 
64757a83faaSBarry Smith   ierr = PetscLogEventBegin(SNESLINESEARCH_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
648bf7f4e0aSPeter Brune 
649bf7f4e0aSPeter Brune   ierr = (*linesearch->ops->apply)(linesearch);CHKERRQ(ierr);
650bf7f4e0aSPeter Brune 
65157a83faaSBarry Smith   ierr = PetscLogEventEnd(SNESLINESEARCH_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
652bf7f4e0aSPeter Brune 
653f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
654bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
655bf7f4e0aSPeter Brune }
656bf7f4e0aSPeter Brune 
657f40b411bSPeter Brune /*@
658f1c6b773SPeter Brune    SNESLineSearchDestroy - Destroys the line search instance.
659f40b411bSPeter Brune 
660f1c6b773SPeter Brune    Collective on SNESLineSearch
661f40b411bSPeter Brune 
662f40b411bSPeter Brune    Input Parameters:
6638e557f58SPeter Brune .  linesearch - The linesearch context
664f40b411bSPeter Brune 
66584238204SBarry Smith    Level: developer
666f40b411bSPeter Brune 
66778bcb3b5SPeter Brune .keywords: SNESLineSearch, Destroy
668f40b411bSPeter Brune 
66984238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchCreate(), SNESLineSearchReset(), SNESDestroy()
670f40b411bSPeter Brune @*/
671bf388a1fSBarry Smith PetscErrorCode SNESLineSearchDestroy(SNESLineSearch * linesearch)
672bf388a1fSBarry Smith {
673bf7f4e0aSPeter Brune   PetscErrorCode ierr;
674bf388a1fSBarry Smith 
675bf7f4e0aSPeter Brune   PetscFunctionBegin;
676bf7f4e0aSPeter Brune   if (!*linesearch) PetscFunctionReturn(0);
677f1c6b773SPeter Brune   PetscValidHeaderSpecific((*linesearch),SNESLINESEARCH_CLASSID,1);
678bf7f4e0aSPeter Brune   if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = 0; PetscFunctionReturn(0);}
679e04113cfSBarry Smith   ierr = PetscObjectSAWsViewOff((PetscObject)*linesearch);CHKERRQ(ierr);
68022d28d08SBarry Smith   ierr = SNESLineSearchReset(*linesearch);CHKERRQ(ierr);
681f5af7f23SKarl Rupp   if ((*linesearch)->ops->destroy) (*linesearch)->ops->destroy(*linesearch);
682bf7f4e0aSPeter Brune   ierr = PetscViewerDestroy(&(*linesearch)->monitor);CHKERRQ(ierr);
683dcf2fd19SBarry Smith   ierr = SNESLineSearchMonitorCancel((*linesearch));CHKERRQ(ierr);
684e7058c64SPeter Brune   ierr = PetscHeaderDestroy(linesearch);CHKERRQ(ierr);
685bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
686bf7f4e0aSPeter Brune }
687bf7f4e0aSPeter Brune 
688f40b411bSPeter Brune /*@
689dcf2fd19SBarry Smith    SNESLineSearchSetDefaultMonitor - Turns on/off printing useful information and debugging output about the line search.
690bf7f4e0aSPeter Brune 
691bf7f4e0aSPeter Brune    Input Parameters:
692dcf2fd19SBarry Smith +  linesearch - the linesearch object
693dcf2fd19SBarry Smith -  viewer - an ASCII PetscViewer or NULL to turn off monitor
694bf7f4e0aSPeter Brune 
695dcf2fd19SBarry Smith    Logically Collective on SNESLineSearch
696bf7f4e0aSPeter Brune 
697bf7f4e0aSPeter Brune    Options Database:
698dcf2fd19SBarry Smith .   -snes_linesearch_monitor [:filename] - enables the monitor
699bf7f4e0aSPeter Brune 
700bf7f4e0aSPeter Brune    Level: intermediate
701bf7f4e0aSPeter Brune 
702dcf2fd19SBarry Smith    Developer Note: This monitor is implemented differently than the other SNESLineSearchMonitors that are set with
703d12e167eSBarry Smith      SNESLineSearchMonitorSet() since it is called in many locations of the line search routines to display aspects of the
704d12e167eSBarry Smith      line search that are not visible to the other monitors.
705dcf2fd19SBarry Smith 
70684238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchGetDefaultMonitor(), PetscViewer, SNESLineSearchSetMonitor()
707bf7f4e0aSPeter Brune @*/
708dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchSetDefaultMonitor(SNESLineSearch linesearch, PetscViewer viewer)
709bf7f4e0aSPeter Brune {
710bf7f4e0aSPeter Brune   PetscErrorCode ierr;
711bf388a1fSBarry Smith 
712bf7f4e0aSPeter Brune   PetscFunctionBegin;
713dcf2fd19SBarry Smith   if (viewer) {ierr = PetscObjectReference((PetscObject)viewer);CHKERRQ(ierr);}
714bf7f4e0aSPeter Brune   ierr = PetscViewerDestroy(&linesearch->monitor);CHKERRQ(ierr);
715dcf2fd19SBarry Smith   linesearch->monitor = viewer;
716bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
717bf7f4e0aSPeter Brune }
718bf7f4e0aSPeter Brune 
719f40b411bSPeter Brune /*@
720dcf2fd19SBarry Smith    SNESLineSearchGetDefaultMonitor - Gets the PetscViewer instance for the line search monitor.
7216a388c36SPeter Brune 
722f190f2fcSBarry Smith    Input Parameter:
7238e557f58SPeter Brune .  linesearch - linesearch context
724f40b411bSPeter Brune 
725f190f2fcSBarry Smith    Output Parameter:
7268e557f58SPeter Brune .  monitor - monitor context
727f40b411bSPeter Brune 
728f40b411bSPeter Brune    Logically Collective on SNES
729f40b411bSPeter Brune 
730f40b411bSPeter Brune    Options Database Keys:
7318e557f58SPeter Brune .   -snes_linesearch_monitor - enables the monitor
732f40b411bSPeter Brune 
733f40b411bSPeter Brune    Level: intermediate
734f40b411bSPeter Brune 
73584238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchSetDefaultMonitor(), PetscViewer
736f40b411bSPeter Brune @*/
737dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchGetDefaultMonitor(SNESLineSearch linesearch, PetscViewer *monitor)
7386a388c36SPeter Brune {
7396a388c36SPeter Brune   PetscFunctionBegin;
740f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
7416a388c36SPeter Brune   if (monitor) {
7426a388c36SPeter Brune     PetscValidPointer(monitor, 2);
7436a388c36SPeter Brune     *monitor = linesearch->monitor;
7446a388c36SPeter Brune   }
7456a388c36SPeter Brune   PetscFunctionReturn(0);
7466a388c36SPeter Brune }
7476a388c36SPeter Brune 
748dcf2fd19SBarry Smith /*@C
749dcf2fd19SBarry Smith    SNESLineSearchMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
750dcf2fd19SBarry Smith 
751dcf2fd19SBarry Smith    Collective on SNESLineSearch
752dcf2fd19SBarry Smith 
753dcf2fd19SBarry Smith    Input Parameters:
754dcf2fd19SBarry Smith +  ls - LineSearch object you wish to monitor
755dcf2fd19SBarry Smith .  name - the monitor type one is seeking
756dcf2fd19SBarry Smith .  help - message indicating what monitoring is done
757dcf2fd19SBarry Smith .  manual - manual page for the monitor
758dcf2fd19SBarry Smith .  monitor - the monitor function
759dcf2fd19SBarry 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
760dcf2fd19SBarry Smith 
761dcf2fd19SBarry Smith    Level: developer
762dcf2fd19SBarry Smith 
763dcf2fd19SBarry Smith .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
764dcf2fd19SBarry Smith           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
765dcf2fd19SBarry Smith           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
766dcf2fd19SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
767dcf2fd19SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
768dcf2fd19SBarry Smith           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
769dcf2fd19SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
770dcf2fd19SBarry Smith @*/
771d12e167eSBarry Smith PetscErrorCode  SNESLineSearchMonitorSetFromOptions(SNESLineSearch ls,const char name[],const char help[], const char manual[],PetscErrorCode (*monitor)(SNESLineSearch,PetscViewerAndFormat*),PetscErrorCode (*monitorsetup)(SNESLineSearch,PetscViewerAndFormat*))
772dcf2fd19SBarry Smith {
773dcf2fd19SBarry Smith   PetscErrorCode    ierr;
774dcf2fd19SBarry Smith   PetscViewer       viewer;
775dcf2fd19SBarry Smith   PetscViewerFormat format;
776dcf2fd19SBarry Smith   PetscBool         flg;
777dcf2fd19SBarry Smith 
778dcf2fd19SBarry Smith   PetscFunctionBegin;
77916413a6aSBarry Smith   ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)ls),((PetscObject) ls)->options,((PetscObject)ls)->prefix,name,&viewer,&format,&flg);CHKERRQ(ierr);
780dcf2fd19SBarry Smith   if (flg) {
781d12e167eSBarry Smith     PetscViewerAndFormat *vf;
782d12e167eSBarry Smith     ierr = PetscViewerAndFormatCreate(viewer,format,&vf);CHKERRQ(ierr);
783d12e167eSBarry Smith     ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr);
784dcf2fd19SBarry Smith     if (monitorsetup) {
785d12e167eSBarry Smith       ierr = (*monitorsetup)(ls,vf);CHKERRQ(ierr);
786dcf2fd19SBarry Smith     }
787d12e167eSBarry Smith     ierr = SNESLineSearchMonitorSet(ls,(PetscErrorCode (*)(SNESLineSearch,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);CHKERRQ(ierr);
788dcf2fd19SBarry Smith   }
789dcf2fd19SBarry Smith   PetscFunctionReturn(0);
790dcf2fd19SBarry Smith }
791dcf2fd19SBarry Smith 
792f40b411bSPeter Brune /*@
793f1c6b773SPeter Brune    SNESLineSearchSetFromOptions - Sets options for the line search
794f40b411bSPeter Brune 
795f40b411bSPeter Brune    Input Parameters:
7968e557f58SPeter Brune .  linesearch - linesearch context
797f40b411bSPeter Brune 
798f40b411bSPeter Brune    Options Database Keys:
799d4c6564cSPatrick Farrell + -snes_linesearch_type <type> - basic, bt, l2, cp, nleqerr, shell
8003c7d6663SPeter Brune . -snes_linesearch_order <order> - 1, 2, 3.  Most types only support certain orders (bt supports 2 or 3)
8011fe24845SBarry Smith . -snes_linesearch_norms   - Turn on/off the linesearch norms for the basic linesearch typem (SNESLineSearchSetComputeNorms())
80271eef1aeSPeter Brune . -snes_linesearch_minlambda - The minimum step length
8031a9b3a06SPeter Brune . -snes_linesearch_maxstep - The maximum step size
8041a9b3a06SPeter Brune . -snes_linesearch_rtol - Relative tolerance for iterative line searches
8051a9b3a06SPeter Brune . -snes_linesearch_atol - Absolute tolerance for iterative line searches
8061a9b3a06SPeter Brune . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches
8071a9b3a06SPeter Brune . -snes_linesearch_max_it - The number of iterations for iterative line searches
808dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches
809dcf2fd19SBarry Smith . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine
8108e557f58SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter
811cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess.
8121a9b3a06SPeter Brune . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method
8131a9b3a06SPeter Brune - -snes_linesearch_precheck_picard_angle - Angle used in picard precheck method
814f40b411bSPeter Brune 
815f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
816f40b411bSPeter Brune 
817f40b411bSPeter Brune    Level: intermediate
818f40b411bSPeter Brune 
8191fe24845SBarry Smith .seealso: SNESLineSearchCreate(), SNESLineSearchSetOrder(), SNESLineSearchSetType(), SNESLineSearchSetTolerances(), SNESLineSearchSetDamping(), SNESLineSearchPreCheckPicard(),
8201fe24845SBarry Smith           SNESLineSearchType, SNESLineSearchSetComputeNorms()
821f40b411bSPeter Brune @*/
822bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch)
823bf388a1fSBarry Smith {
824bf7f4e0aSPeter Brune   PetscErrorCode    ierr;
8251a4f838cSPeter Brune   const char        *deft = SNESLINESEARCHBASIC;
826bf7f4e0aSPeter Brune   char              type[256];
827bf7f4e0aSPeter Brune   PetscBool         flg, set;
828dcf2fd19SBarry Smith   PetscViewer       viewer;
829bf388a1fSBarry Smith 
830bf7f4e0aSPeter Brune   PetscFunctionBegin;
8310f51fdf8SToby Isaac   ierr = SNESLineSearchRegisterAll();CHKERRQ(ierr);
832bf7f4e0aSPeter Brune 
833bf7f4e0aSPeter Brune   ierr = PetscObjectOptionsBegin((PetscObject)linesearch);CHKERRQ(ierr);
834f5af7f23SKarl Rupp   if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name;
835a264d7a6SBarry Smith   ierr = PetscOptionsFList("-snes_linesearch_type","Linesearch type","SNESLineSearchSetType",SNESLineSearchList,deft,type,256,&flg);CHKERRQ(ierr);
836bf7f4e0aSPeter Brune   if (flg) {
837f1c6b773SPeter Brune     ierr = SNESLineSearchSetType(linesearch,type);CHKERRQ(ierr);
838bf7f4e0aSPeter Brune   } else if (!((PetscObject)linesearch)->type_name) {
839f1c6b773SPeter Brune     ierr = SNESLineSearchSetType(linesearch,deft);CHKERRQ(ierr);
840bf7f4e0aSPeter Brune   }
841bf7f4e0aSPeter Brune 
84216413a6aSBarry Smith   ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)linesearch),((PetscObject) linesearch)->options,((PetscObject)linesearch)->prefix,"-snes_linesearch_monitor",&viewer,NULL,&set);CHKERRQ(ierr);
843dcf2fd19SBarry Smith   if (set) {
844dcf2fd19SBarry Smith     ierr = SNESLineSearchSetDefaultMonitor(linesearch,viewer);CHKERRQ(ierr);
845dcf2fd19SBarry Smith     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
846dcf2fd19SBarry Smith   }
847dcf2fd19SBarry Smith   ierr = SNESLineSearchMonitorSetFromOptions(linesearch,"-snes_linesearch_monitor_solution_update","View correction at each iteration","SNESLineSearchMonitorSolutionUpdate",SNESLineSearchMonitorSolutionUpdate,NULL);CHKERRQ(ierr);
848bf7f4e0aSPeter Brune 
8491a9b3a06SPeter Brune   /* tolerances */
85094ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_minlambda","Minimum step length","SNESLineSearchSetTolerances",linesearch->steptol,&linesearch->steptol,NULL);CHKERRQ(ierr);
85194ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_maxstep","Maximum step size","SNESLineSearchSetTolerances",linesearch->maxstep,&linesearch->maxstep,NULL);CHKERRQ(ierr);
85294ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_rtol","Relative tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->rtol,&linesearch->rtol,NULL);CHKERRQ(ierr);
85394ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_atol","Absolute tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->atol,&linesearch->atol,NULL);CHKERRQ(ierr);
85494ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_ltol","Change in lambda tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->ltol,&linesearch->ltol,NULL);CHKERRQ(ierr);
85594ae4db5SBarry Smith   ierr = PetscOptionsInt("-snes_linesearch_max_it","Maximum iterations for iterative line searches","SNESLineSearchSetTolerances",linesearch->max_its,&linesearch->max_its,NULL);CHKERRQ(ierr);
8567a35526eSPeter Brune 
8571a9b3a06SPeter Brune   /* damping parameters */
85894ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_damping","Line search damping and initial step guess","SNESLineSearchSetDamping",linesearch->damping,&linesearch->damping,NULL);CHKERRQ(ierr);
8591a9b3a06SPeter Brune 
86094ae4db5SBarry Smith   ierr = PetscOptionsBool("-snes_linesearch_keeplambda","Use previous lambda as damping","SNESLineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,NULL);CHKERRQ(ierr);
8611a9b3a06SPeter Brune 
8621a9b3a06SPeter Brune   /* precheck */
8637a35526eSPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_precheck_picard","Use a correction that sometimes improves convergence of Picard iteration","SNESLineSearchPreCheckPicard",flg,&flg,&set);CHKERRQ(ierr);
8647a35526eSPeter Brune   if (set) {
8657a35526eSPeter Brune     if (flg) {
8667a35526eSPeter Brune       linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */
867f5af7f23SKarl Rupp 
8687a35526eSPeter Brune       ierr = PetscOptionsReal("-snes_linesearch_precheck_picard_angle","Maximum angle at which to activate the correction",
8690298fd71SBarry Smith                               "none",linesearch->precheck_picard_angle,&linesearch->precheck_picard_angle,NULL);CHKERRQ(ierr);
8707a35526eSPeter Brune       ierr = SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&linesearch->precheck_picard_angle);CHKERRQ(ierr);
8717a35526eSPeter Brune     } else {
8720298fd71SBarry Smith       ierr = SNESLineSearchSetPreCheck(linesearch,NULL,NULL);CHKERRQ(ierr);
8737a35526eSPeter Brune     }
8747a35526eSPeter Brune   }
87594ae4db5SBarry Smith   ierr = PetscOptionsInt("-snes_linesearch_order","Order of approximation used in the line search","SNESLineSearchSetOrder",linesearch->order,&linesearch->order,NULL);CHKERRQ(ierr);
87694ae4db5SBarry Smith   ierr = PetscOptionsBool("-snes_linesearch_norms","Compute final norms in line search","SNESLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,NULL);CHKERRQ(ierr);
8777a35526eSPeter Brune 
8785a9b6599SPeter Brune   if (linesearch->ops->setfromoptions) {
879e55864a3SBarry Smith     (*linesearch->ops->setfromoptions)(PetscOptionsObject,linesearch);CHKERRQ(ierr);
8805a9b6599SPeter Brune   }
8815a9b6599SPeter Brune 
8820633abcbSJed Brown   ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)linesearch);CHKERRQ(ierr);
883bf7f4e0aSPeter Brune   ierr = PetscOptionsEnd();CHKERRQ(ierr);
884bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
885bf7f4e0aSPeter Brune }
886bf7f4e0aSPeter Brune 
887f40b411bSPeter Brune /*@
888f190f2fcSBarry Smith    SNESLineSearchView - Prints useful information about the line search
889f40b411bSPeter Brune 
890f40b411bSPeter Brune    Input Parameters:
8918e557f58SPeter Brune .  linesearch - linesearch context
892f40b411bSPeter Brune 
893f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
894f40b411bSPeter Brune 
895f40b411bSPeter Brune    Level: intermediate
896f40b411bSPeter Brune 
897dcf2fd19SBarry Smith .seealso: SNESLineSearchCreate()
898f40b411bSPeter Brune @*/
899bf388a1fSBarry Smith PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer)
900bf388a1fSBarry Smith {
9017f1410a3SPeter Brune   PetscErrorCode ierr;
9027f1410a3SPeter Brune   PetscBool      iascii;
903bf388a1fSBarry Smith 
904bf7f4e0aSPeter Brune   PetscFunctionBegin;
9057f1410a3SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
9067f1410a3SPeter Brune   if (!viewer) {
907ce94432eSBarry Smith     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch),&viewer);CHKERRQ(ierr);
9087f1410a3SPeter Brune   }
9097f1410a3SPeter Brune   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
9107f1410a3SPeter Brune   PetscCheckSameComm(linesearch,1,viewer,2);
911f40b411bSPeter Brune 
912251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
9137f1410a3SPeter Brune   if (iascii) {
914dae58748SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)linesearch,viewer);CHKERRQ(ierr);
9157f1410a3SPeter Brune     if (linesearch->ops->view) {
9167f1410a3SPeter Brune       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
9177f1410a3SPeter Brune       ierr = (*linesearch->ops->view)(linesearch,viewer);CHKERRQ(ierr);
9187f1410a3SPeter Brune       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
9197f1410a3SPeter Brune     }
920c69d1a72SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep,(double)linesearch->steptol);CHKERRQ(ierr);
921c69d1a72SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol,(double)linesearch->atol,(double)linesearch->ltol);CHKERRQ(ierr);
9227f1410a3SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D\n", linesearch->max_its);CHKERRQ(ierr);
9236b2b7091SBarry Smith     if (linesearch->ops->precheck) {
9246b2b7091SBarry Smith       if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) {
9257f1410a3SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  using precheck step to speed up Picard convergence\n", linesearch->max_its);CHKERRQ(ierr);
9267f1410a3SPeter Brune       } else {
9277f1410a3SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  using user-defined precheck step\n", linesearch->max_its);CHKERRQ(ierr);
9287f1410a3SPeter Brune       }
9297f1410a3SPeter Brune     }
9306b2b7091SBarry Smith     if (linesearch->ops->postcheck) {
9317f1410a3SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer,"  using user-defined postcheck step\n", linesearch->max_its);CHKERRQ(ierr);
9327f1410a3SPeter Brune     }
9337f1410a3SPeter Brune   }
934bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
935bf7f4e0aSPeter Brune }
936bf7f4e0aSPeter Brune 
937ea5d4fccSPeter Brune /*@C
938f1c6b773SPeter Brune    SNESLineSearchSetType - Sets the linesearch type
939f40b411bSPeter Brune 
940f190f2fcSBarry Smith    Logically Collective on SNESLineSearch
941f190f2fcSBarry Smith 
942f40b411bSPeter Brune    Input Parameters:
9438e557f58SPeter Brune +  linesearch - linesearch context
944f40b411bSPeter Brune -  type - The type of line search to be used
945f40b411bSPeter Brune 
946cd7522eaSPeter Brune    Available Types:
9471fe24845SBarry Smith +  SNESLINESEARCHBASIC - Simple damping line search, defaults to using the full Newton step
9481fe24845SBarry Smith .  SNESLINESEARCHBT - Backtracking line search over the L2 norm of the function
9491fe24845SBarry Smith .  SNESLINESEARCHL2 - Secant line search over the L2 norm of the function
9501fe24845SBarry Smith .  SNESLINESEARCHCP - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x)
9511fe24845SBarry Smith .  SNESLINESEARCHNLEQERR - Affine-covariant error-oriented linesearch
9521fe24845SBarry Smith -  SNESLINESEARCHSHELL - User provided SNESLineSearch implementation
9531fe24845SBarry Smith 
9541fe24845SBarry Smith    Options Database:
9551fe24845SBarry Smith .  -snes_linesearch_type <type> - basic, bt, l2, cp, nleqerr, shell
956cd7522eaSPeter Brune 
957f40b411bSPeter Brune    Level: intermediate
958f40b411bSPeter Brune 
9591fe24845SBarry Smith .seealso: SNESLineSearchCreate(), SNESLineSearchType, SNESLineSearchSetFromOptions()
960f40b411bSPeter Brune @*/
96119fd82e9SBarry Smith PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type)
962bf7f4e0aSPeter Brune {
963f1c6b773SPeter Brune   PetscErrorCode ierr,(*r)(SNESLineSearch);
964bf7f4e0aSPeter Brune   PetscBool      match;
965bf7f4e0aSPeter Brune 
966bf7f4e0aSPeter Brune   PetscFunctionBegin;
967f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
968bf7f4e0aSPeter Brune   PetscValidCharPointer(type,2);
969bf7f4e0aSPeter Brune 
970251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)linesearch,type,&match);CHKERRQ(ierr);
971bf7f4e0aSPeter Brune   if (match) PetscFunctionReturn(0);
972bf7f4e0aSPeter Brune 
9731c9cd337SJed Brown   ierr = PetscFunctionListFind(SNESLineSearchList,type,&r);CHKERRQ(ierr);
974bf7f4e0aSPeter Brune   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type);
975bf7f4e0aSPeter Brune   /* Destroy the previous private linesearch context */
976bf7f4e0aSPeter Brune   if (linesearch->ops->destroy) {
977bf7f4e0aSPeter Brune     ierr = (*(linesearch)->ops->destroy)(linesearch);CHKERRQ(ierr);
978f5af7f23SKarl Rupp 
9790298fd71SBarry Smith     linesearch->ops->destroy = NULL;
980bf7f4e0aSPeter Brune   }
981f1c6b773SPeter Brune   /* Reinitialize function pointers in SNESLineSearchOps structure */
982bf7f4e0aSPeter Brune   linesearch->ops->apply          = 0;
983bf7f4e0aSPeter Brune   linesearch->ops->view           = 0;
984bf7f4e0aSPeter Brune   linesearch->ops->setfromoptions = 0;
985bf7f4e0aSPeter Brune   linesearch->ops->destroy        = 0;
986bf7f4e0aSPeter Brune 
987bf7f4e0aSPeter Brune   ierr = PetscObjectChangeTypeName((PetscObject)linesearch,type);CHKERRQ(ierr);
988bf7f4e0aSPeter Brune   ierr = (*r)(linesearch);CHKERRQ(ierr);
989bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
990bf7f4e0aSPeter Brune }
991bf7f4e0aSPeter Brune 
992f40b411bSPeter Brune /*@
99378bcb3b5SPeter Brune    SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation.
994f40b411bSPeter Brune 
995f40b411bSPeter Brune    Input Parameters:
9968e557f58SPeter Brune +  linesearch - linesearch context
997f40b411bSPeter Brune -  snes - The snes instance
998f40b411bSPeter Brune 
99978bcb3b5SPeter Brune    Level: developer
100078bcb3b5SPeter Brune 
100178bcb3b5SPeter Brune    Notes:
1002f190f2fcSBarry Smith    This happens automatically when the line search is obtained/created with
10037601faf0SJed Brown    SNESGetLineSearch().  This routine is therefore mainly called within SNES
100478bcb3b5SPeter Brune    implementations.
1005f40b411bSPeter Brune 
10068141a3b9SPeter Brune    Level: developer
1007f40b411bSPeter Brune 
1008cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
1009f40b411bSPeter Brune @*/
1010bf388a1fSBarry Smith PetscErrorCode  SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes)
1011bf388a1fSBarry Smith {
1012bf7f4e0aSPeter Brune   PetscFunctionBegin;
1013f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1014bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,2);
1015bf7f4e0aSPeter Brune   linesearch->snes = snes;
1016bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1017bf7f4e0aSPeter Brune }
1018bf7f4e0aSPeter Brune 
1019f40b411bSPeter Brune /*@
10208141a3b9SPeter Brune    SNESLineSearchGetSNES - Gets the SNES instance associated with the line search.
10218141a3b9SPeter Brune    Having an associated SNES is necessary because most line search implementations must be able to
10228141a3b9SPeter Brune    evaluate the function using SNESComputeFunction() for the associated SNES.  This routine
10238141a3b9SPeter Brune    is used in the line search implementations when one must get this associated SNES instance.
1024f40b411bSPeter Brune 
1025f40b411bSPeter Brune    Input Parameters:
10268e557f58SPeter Brune .  linesearch - linesearch context
1027f40b411bSPeter Brune 
1028f40b411bSPeter Brune    Output Parameters:
1029f40b411bSPeter Brune .  snes - The snes instance
1030f40b411bSPeter Brune 
10318141a3b9SPeter Brune    Level: developer
1032f40b411bSPeter Brune 
1033cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
1034f40b411bSPeter Brune @*/
1035bf388a1fSBarry Smith PetscErrorCode  SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes)
1036bf388a1fSBarry Smith {
1037bf7f4e0aSPeter Brune   PetscFunctionBegin;
1038f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
10396a388c36SPeter Brune   PetscValidPointer(snes, 2);
1040bf7f4e0aSPeter Brune   *snes = linesearch->snes;
1041bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1042bf7f4e0aSPeter Brune }
1043bf7f4e0aSPeter Brune 
1044f40b411bSPeter Brune /*@
1045f1c6b773SPeter Brune    SNESLineSearchGetLambda - Gets the last linesearch steplength discovered.
1046f40b411bSPeter Brune 
1047f40b411bSPeter Brune    Input Parameters:
10488e557f58SPeter Brune .  linesearch - linesearch context
1049f40b411bSPeter Brune 
1050f40b411bSPeter Brune    Output Parameters:
1051cd7522eaSPeter Brune .  lambda - The last steplength computed during SNESLineSearchApply()
1052f40b411bSPeter Brune 
105378bcb3b5SPeter Brune    Level: advanced
105478bcb3b5SPeter Brune 
10558e557f58SPeter Brune    Notes:
10568e557f58SPeter Brune    This is useful in methods where the solver is ill-scaled and
105778bcb3b5SPeter Brune    requires some adaptive notion of the difference in scale between the
105878bcb3b5SPeter Brune    solution and the function.  For instance, SNESQN may be scaled by the
105978bcb3b5SPeter Brune    line search lambda using the argument -snes_qn_scaling ls.
106078bcb3b5SPeter Brune 
1061cd7522eaSPeter Brune .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetDamping(), SNESLineSearchApply()
1062f40b411bSPeter Brune @*/
1063f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda)
10646a388c36SPeter Brune {
10656a388c36SPeter Brune   PetscFunctionBegin;
1066f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
10676a388c36SPeter Brune   PetscValidPointer(lambda, 2);
10686a388c36SPeter Brune   *lambda = linesearch->lambda;
10696a388c36SPeter Brune   PetscFunctionReturn(0);
10706a388c36SPeter Brune }
10716a388c36SPeter Brune 
1072f40b411bSPeter Brune /*@
1073f1c6b773SPeter Brune    SNESLineSearchSetLambda - Sets the linesearch steplength.
1074f40b411bSPeter Brune 
1075f40b411bSPeter Brune    Input Parameters:
10768e557f58SPeter Brune +  linesearch - linesearch context
1077f40b411bSPeter Brune -  lambda - The last steplength.
1078f40b411bSPeter Brune 
1079cd7522eaSPeter Brune    Notes:
1080f190f2fcSBarry Smith    This routine is typically used within implementations of SNESLineSearchApply()
1081cd7522eaSPeter Brune    to set the final steplength.  This routine (and SNESLineSearchGetLambda()) were
1082cd7522eaSPeter Brune    added in order to facilitate Quasi-Newton methods that use the previous steplength
1083cd7522eaSPeter Brune    as an inner scaling parameter.
1084cd7522eaSPeter Brune 
108578bcb3b5SPeter Brune    Level: advanced
1086f40b411bSPeter Brune 
1087f1c6b773SPeter Brune .seealso: SNESLineSearchGetLambda()
1088f40b411bSPeter Brune @*/
1089f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
10906a388c36SPeter Brune {
10916a388c36SPeter Brune   PetscFunctionBegin;
1092f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
10936a388c36SPeter Brune   linesearch->lambda = lambda;
10946a388c36SPeter Brune   PetscFunctionReturn(0);
10956a388c36SPeter Brune }
10966a388c36SPeter Brune 
1097f40b411bSPeter Brune /*@
10983c7d6663SPeter Brune    SNESLineSearchGetTolerances - Gets the tolerances for the linesearch.  These include
109978bcb3b5SPeter Brune    tolerances for the relative and absolute change in the function norm, the change
110078bcb3b5SPeter Brune    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
110178bcb3b5SPeter Brune    and the maximum number of iterations the line search procedure may take.
1102f40b411bSPeter Brune 
1103f40b411bSPeter Brune    Input Parameters:
11048e557f58SPeter Brune .  linesearch - linesearch context
1105f40b411bSPeter Brune 
1106f40b411bSPeter Brune    Output Parameters:
1107516fe3c3SPeter Brune +  steptol - The minimum steplength
11086cc8e53bSPeter Brune .  maxstep - The maximum steplength
1109516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
1110516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
1111516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
1112516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
1113f40b411bSPeter Brune 
111478bcb3b5SPeter Brune    Level: intermediate
111578bcb3b5SPeter Brune 
111678bcb3b5SPeter Brune    Notes:
111778bcb3b5SPeter Brune    Different line searches may implement these parameters slightly differently as
11183c7d6663SPeter Brune    the type requires.
1119516fe3c3SPeter Brune 
1120f1c6b773SPeter Brune .seealso: SNESLineSearchSetTolerances()
1121f40b411bSPeter Brune @*/
1122f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
11236a388c36SPeter Brune {
11246a388c36SPeter Brune   PetscFunctionBegin;
1125f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1126516fe3c3SPeter Brune   if (steptol) {
11276a388c36SPeter Brune     PetscValidPointer(steptol, 2);
11286a388c36SPeter Brune     *steptol = linesearch->steptol;
1129516fe3c3SPeter Brune   }
1130516fe3c3SPeter Brune   if (maxstep) {
1131516fe3c3SPeter Brune     PetscValidPointer(maxstep, 3);
1132516fe3c3SPeter Brune     *maxstep = linesearch->maxstep;
1133516fe3c3SPeter Brune   }
1134516fe3c3SPeter Brune   if (rtol) {
1135516fe3c3SPeter Brune     PetscValidPointer(rtol, 4);
1136516fe3c3SPeter Brune     *rtol = linesearch->rtol;
1137516fe3c3SPeter Brune   }
1138516fe3c3SPeter Brune   if (atol) {
1139516fe3c3SPeter Brune     PetscValidPointer(atol, 5);
1140516fe3c3SPeter Brune     *atol = linesearch->atol;
1141516fe3c3SPeter Brune   }
1142516fe3c3SPeter Brune   if (ltol) {
1143516fe3c3SPeter Brune     PetscValidPointer(ltol, 6);
1144516fe3c3SPeter Brune     *ltol = linesearch->ltol;
1145516fe3c3SPeter Brune   }
1146516fe3c3SPeter Brune   if (max_its) {
1147516fe3c3SPeter Brune     PetscValidPointer(max_its, 7);
1148516fe3c3SPeter Brune     *max_its = linesearch->max_its;
1149516fe3c3SPeter Brune   }
11506a388c36SPeter Brune   PetscFunctionReturn(0);
11516a388c36SPeter Brune }
11526a388c36SPeter Brune 
1153f40b411bSPeter Brune /*@
11543c7d6663SPeter Brune    SNESLineSearchSetTolerances -  Gets the tolerances for the linesearch.  These include
115578bcb3b5SPeter Brune    tolerances for the relative and absolute change in the function norm, the change
115678bcb3b5SPeter Brune    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
115778bcb3b5SPeter Brune    and the maximum number of iterations the line search procedure may take.
1158f40b411bSPeter Brune 
1159f40b411bSPeter Brune    Input Parameters:
11608e557f58SPeter Brune +  linesearch - linesearch context
1161516fe3c3SPeter Brune .  steptol - The minimum steplength
11626cc8e53bSPeter Brune .  maxstep - The maximum steplength
1163516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
1164516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
1165516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
1166516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
1167f40b411bSPeter Brune 
116878bcb3b5SPeter Brune    Notes:
11693c7d6663SPeter Brune    The user may choose to not set any of the tolerances using PETSC_DEFAULT in
117078bcb3b5SPeter Brune    place of an argument.
1171f40b411bSPeter Brune 
117278bcb3b5SPeter Brune    Level: intermediate
1173516fe3c3SPeter Brune 
1174f1c6b773SPeter Brune .seealso: SNESLineSearchGetTolerances()
1175f40b411bSPeter Brune @*/
1176f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its)
11776a388c36SPeter Brune {
11786a388c36SPeter Brune   PetscFunctionBegin;
1179f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1180d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,steptol,2);
1181d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,maxstep,3);
1182d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,rtol,4);
1183d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,atol,5);
1184d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,ltol,6);
1185d3952184SSatish Balay   PetscValidLogicalCollectiveInt(linesearch,max_its,7);
1186d3952184SSatish Balay 
1187d3952184SSatish Balay   if (steptol!= PETSC_DEFAULT) {
1188ce94432eSBarry Smith     if (steptol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Minimum step length %14.12e must be non-negative",(double)steptol);
11896a388c36SPeter Brune     linesearch->steptol = steptol;
1190d3952184SSatish Balay   }
1191d3952184SSatish Balay 
1192d3952184SSatish Balay   if (maxstep!= PETSC_DEFAULT) {
1193ce94432eSBarry Smith     if (maxstep < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum step length %14.12e must be non-negative",(double)maxstep);
1194516fe3c3SPeter Brune     linesearch->maxstep = maxstep;
1195d3952184SSatish Balay   }
1196d3952184SSatish Balay 
1197d3952184SSatish Balay   if (rtol != PETSC_DEFAULT) {
1198ce94432eSBarry Smith     if (rtol < 0.0 || 1.0 <= rtol) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Relative tolerance %14.12e must be non-negative and less than 1.0",(double)rtol);
1199516fe3c3SPeter Brune     linesearch->rtol = rtol;
1200d3952184SSatish Balay   }
1201d3952184SSatish Balay 
1202d3952184SSatish Balay   if (atol != PETSC_DEFAULT) {
1203ce94432eSBarry Smith     if (atol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %14.12e must be non-negative",(double)atol);
1204516fe3c3SPeter Brune     linesearch->atol = atol;
1205d3952184SSatish Balay   }
1206d3952184SSatish Balay 
1207d3952184SSatish Balay   if (ltol != PETSC_DEFAULT) {
1208ce94432eSBarry Smith     if (ltol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Labmda tolerance %14.12e must be non-negative",(double)ltol);
1209516fe3c3SPeter Brune     linesearch->ltol = ltol;
1210d3952184SSatish Balay   }
1211d3952184SSatish Balay 
1212d3952184SSatish Balay   if (max_its != PETSC_DEFAULT) {
1213ce94432eSBarry Smith     if (max_its < 0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",max_its);
1214516fe3c3SPeter Brune     linesearch->max_its = max_its;
1215d3952184SSatish Balay   }
12166a388c36SPeter Brune   PetscFunctionReturn(0);
12176a388c36SPeter Brune }
12186a388c36SPeter Brune 
1219f40b411bSPeter Brune /*@
1220f1c6b773SPeter Brune    SNESLineSearchGetDamping - Gets the line search damping parameter.
1221f40b411bSPeter Brune 
1222f40b411bSPeter Brune    Input Parameters:
12238e557f58SPeter Brune .  linesearch - linesearch context
1224f40b411bSPeter Brune 
1225f40b411bSPeter Brune    Output Parameters:
12268e557f58SPeter Brune .  damping - The damping parameter
1227f40b411bSPeter Brune 
122878bcb3b5SPeter Brune    Level: advanced
1229f40b411bSPeter Brune 
123078bcb3b5SPeter Brune .seealso: SNESLineSearchGetStepTolerance(), SNESQN
1231f40b411bSPeter Brune @*/
1232f40b411bSPeter Brune 
1233f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping)
12346a388c36SPeter Brune {
12356a388c36SPeter Brune   PetscFunctionBegin;
1236f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
12376a388c36SPeter Brune   PetscValidPointer(damping, 2);
12386a388c36SPeter Brune   *damping = linesearch->damping;
12396a388c36SPeter Brune   PetscFunctionReturn(0);
12406a388c36SPeter Brune }
12416a388c36SPeter Brune 
1242f40b411bSPeter Brune /*@
1243f1c6b773SPeter Brune    SNESLineSearchSetDamping - Sets the line search damping paramter.
1244f40b411bSPeter Brune 
1245f40b411bSPeter Brune    Input Parameters:
124603fd4120SBarry Smith +  linesearch - linesearch context
124703fd4120SBarry Smith -  damping - The damping parameter
1248f40b411bSPeter Brune 
124903fd4120SBarry Smith    Options Database:
125003fd4120SBarry Smith .   -snes_linesearch_damping
1251f40b411bSPeter Brune    Level: intermediate
1252f40b411bSPeter Brune 
1253cd7522eaSPeter Brune    Notes:
1254cd7522eaSPeter Brune    The basic line search merely takes the update step scaled by the damping parameter.
1255cd7522eaSPeter Brune    The use of the damping parameter in the l2 and cp line searches is much more subtle;
125678bcb3b5SPeter Brune    it is used as a starting point in calculating the secant step. However, the eventual
1257cd7522eaSPeter Brune    step may be of greater length than the damping parameter.  In the bt line search it is
1258cd7522eaSPeter Brune    used as the maximum possible step length, as the bt line search only backtracks.
1259cd7522eaSPeter Brune 
1260f1c6b773SPeter Brune .seealso: SNESLineSearchGetDamping()
1261f40b411bSPeter Brune @*/
1262f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping)
12636a388c36SPeter Brune {
12646a388c36SPeter Brune   PetscFunctionBegin;
1265f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
12666a388c36SPeter Brune   linesearch->damping = damping;
12676a388c36SPeter Brune   PetscFunctionReturn(0);
12686a388c36SPeter Brune }
12696a388c36SPeter Brune 
127059405d5eSPeter Brune /*@
127159405d5eSPeter Brune    SNESLineSearchGetOrder - Gets the line search approximation order.
127259405d5eSPeter Brune 
127359405d5eSPeter Brune    Input Parameters:
127478bcb3b5SPeter Brune .  linesearch - linesearch context
127559405d5eSPeter Brune 
127659405d5eSPeter Brune    Output Parameters:
12778e557f58SPeter Brune .  order - The order
127859405d5eSPeter Brune 
127978bcb3b5SPeter Brune    Possible Values for order:
12803c7d6663SPeter Brune +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
12813c7d6663SPeter Brune .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
12823c7d6663SPeter Brune -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
128378bcb3b5SPeter Brune 
128459405d5eSPeter Brune    Level: intermediate
128559405d5eSPeter Brune 
128659405d5eSPeter Brune .seealso: SNESLineSearchSetOrder()
128759405d5eSPeter Brune @*/
128859405d5eSPeter Brune 
1289b000cd8dSPeter Brune PetscErrorCode  SNESLineSearchGetOrder(SNESLineSearch linesearch,PetscInt *order)
129059405d5eSPeter Brune {
129159405d5eSPeter Brune   PetscFunctionBegin;
129259405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
129359405d5eSPeter Brune   PetscValidPointer(order, 2);
129459405d5eSPeter Brune   *order = linesearch->order;
129559405d5eSPeter Brune   PetscFunctionReturn(0);
129659405d5eSPeter Brune }
129759405d5eSPeter Brune 
129859405d5eSPeter Brune /*@
12991f8196a2SJed Brown    SNESLineSearchSetOrder - Sets the maximum order of the polynomial fit used in the line search
130059405d5eSPeter Brune 
130159405d5eSPeter Brune    Input Parameters:
1302*a2b725a8SWilliam Gropp +  linesearch - linesearch context
1303*a2b725a8SWilliam Gropp -  order - The damping parameter
130459405d5eSPeter Brune 
130559405d5eSPeter Brune    Level: intermediate
130659405d5eSPeter Brune 
130778bcb3b5SPeter Brune    Possible Values for order:
13083c7d6663SPeter Brune +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
13093c7d6663SPeter Brune .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
13103c7d6663SPeter Brune -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
131178bcb3b5SPeter Brune 
131259405d5eSPeter Brune    Notes:
131359405d5eSPeter Brune    Variable orders are supported by the following line searches:
131478bcb3b5SPeter Brune +  bt - cubic and quadratic
131578bcb3b5SPeter Brune -  cp - linear and quadratic
131659405d5eSPeter Brune 
13171f8196a2SJed Brown .seealso: SNESLineSearchGetOrder(), SNESLineSearchSetDamping()
131859405d5eSPeter Brune @*/
1319b000cd8dSPeter Brune PetscErrorCode  SNESLineSearchSetOrder(SNESLineSearch linesearch,PetscInt order)
132059405d5eSPeter Brune {
132159405d5eSPeter Brune   PetscFunctionBegin;
132259405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
132359405d5eSPeter Brune   linesearch->order = order;
132459405d5eSPeter Brune   PetscFunctionReturn(0);
132559405d5eSPeter Brune }
132659405d5eSPeter Brune 
1327f40b411bSPeter Brune /*@
1328f1c6b773SPeter Brune    SNESLineSearchGetNorms - Gets the norms for for X, Y, and F.
1329f40b411bSPeter Brune 
1330f40b411bSPeter Brune    Input Parameters:
133178bcb3b5SPeter Brune .  linesearch - linesearch context
1332f40b411bSPeter Brune 
1333f40b411bSPeter Brune    Output Parameters:
1334f40b411bSPeter Brune +  xnorm - The norm of the current solution
1335f40b411bSPeter Brune .  fnorm - The norm of the current function
1336f40b411bSPeter Brune -  ynorm - The norm of the current update
1337f40b411bSPeter Brune 
1338cd7522eaSPeter Brune    Notes:
1339cd7522eaSPeter Brune    This function is mainly called from SNES implementations.
1340cd7522eaSPeter Brune 
134178bcb3b5SPeter Brune    Level: developer
1342f40b411bSPeter Brune 
1343f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms() SNESLineSearchGetVecs()
1344f40b411bSPeter Brune @*/
1345f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm)
1346bf7f4e0aSPeter Brune {
1347bf7f4e0aSPeter Brune   PetscFunctionBegin;
1348f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1349f5af7f23SKarl Rupp   if (xnorm) *xnorm = linesearch->xnorm;
1350f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
1351f5af7f23SKarl Rupp   if (ynorm) *ynorm = linesearch->ynorm;
1352bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1353bf7f4e0aSPeter Brune }
1354bf7f4e0aSPeter Brune 
1355f40b411bSPeter Brune /*@
1356f1c6b773SPeter Brune    SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F.
1357f40b411bSPeter Brune 
1358f40b411bSPeter Brune    Input Parameters:
135978bcb3b5SPeter Brune +  linesearch - linesearch context
1360f40b411bSPeter Brune .  xnorm - The norm of the current solution
1361f40b411bSPeter Brune .  fnorm - The norm of the current function
1362f40b411bSPeter Brune -  ynorm - The norm of the current update
1363f40b411bSPeter Brune 
136478bcb3b5SPeter Brune    Level: advanced
1365f40b411bSPeter Brune 
1366f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1367f40b411bSPeter Brune @*/
1368f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
13696a388c36SPeter Brune {
13706a388c36SPeter Brune   PetscFunctionBegin;
1371f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
13726a388c36SPeter Brune   linesearch->xnorm = xnorm;
13736a388c36SPeter Brune   linesearch->fnorm = fnorm;
13746a388c36SPeter Brune   linesearch->ynorm = ynorm;
13756a388c36SPeter Brune   PetscFunctionReturn(0);
13766a388c36SPeter Brune }
13776a388c36SPeter Brune 
1378f40b411bSPeter Brune /*@
1379f1c6b773SPeter Brune    SNESLineSearchComputeNorms - Computes the norms of X, F, and Y.
1380f40b411bSPeter Brune 
1381f40b411bSPeter Brune    Input Parameters:
138278bcb3b5SPeter Brune .  linesearch - linesearch context
1383f40b411bSPeter Brune 
1384f40b411bSPeter Brune    Options Database Keys:
13858e557f58SPeter Brune .   -snes_linesearch_norms - turn norm computation on or off
1386f40b411bSPeter Brune 
1387f40b411bSPeter Brune    Level: intermediate
1388f40b411bSPeter Brune 
138978bcb3b5SPeter Brune .seealso: SNESLineSearchGetNorms, SNESLineSearchSetNorms(), SNESLineSearchSetComputeNorms()
1390f40b411bSPeter Brune @*/
1391f1c6b773SPeter Brune PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
13926a388c36SPeter Brune {
13936a388c36SPeter Brune   PetscErrorCode ierr;
13949bd66eb0SPeter Brune   SNES           snes;
1395bf388a1fSBarry Smith 
13966a388c36SPeter Brune   PetscFunctionBegin;
13976a388c36SPeter Brune   if (linesearch->norms) {
13989bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
1399f1c6b773SPeter Brune       ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
14009bd66eb0SPeter Brune       ierr = VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
14019bd66eb0SPeter Brune       ierr = VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
14029bd66eb0SPeter Brune       ierr = (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);CHKERRQ(ierr);
14039bd66eb0SPeter Brune     } else {
14046a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_func,   NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
14056a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_sol,    NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
14066a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
14076a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_func,     NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
14086a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_sol,      NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
14096a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_update,   NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
14106a388c36SPeter Brune     }
14119bd66eb0SPeter Brune   }
14126a388c36SPeter Brune   PetscFunctionReturn(0);
14136a388c36SPeter Brune }
14146a388c36SPeter Brune 
14156f263ca3SPeter Brune /*@
14166f263ca3SPeter Brune    SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search.
14176f263ca3SPeter Brune 
14186f263ca3SPeter Brune    Input Parameters:
141978bcb3b5SPeter Brune +  linesearch  - linesearch context
142078bcb3b5SPeter Brune -  flg  - indicates whether or not to compute norms
14216f263ca3SPeter Brune 
14226f263ca3SPeter Brune    Options Database Keys:
14231f8196a2SJed Brown .   -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic linesearch
14246f263ca3SPeter Brune 
14256f263ca3SPeter Brune    Notes:
14261f8196a2SJed Brown    This is most relevant to the SNESLINESEARCHBASIC line search type since most line searches have a stopping criteria involving the norm.
14276f263ca3SPeter Brune 
14286f263ca3SPeter Brune    Level: intermediate
14296f263ca3SPeter Brune 
14301a4f838cSPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetNorms(), SNESLineSearchComputeNorms(), SNESLINESEARCHBASIC
14316f263ca3SPeter Brune @*/
14326f263ca3SPeter Brune PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg)
14336f263ca3SPeter Brune {
14346f263ca3SPeter Brune   PetscFunctionBegin;
14356f263ca3SPeter Brune   linesearch->norms = flg;
14366f263ca3SPeter Brune   PetscFunctionReturn(0);
14376f263ca3SPeter Brune }
14386f263ca3SPeter Brune 
1439f40b411bSPeter Brune /*@
1440f1c6b773SPeter Brune    SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context
1441f40b411bSPeter Brune 
1442f40b411bSPeter Brune    Input Parameters:
144378bcb3b5SPeter Brune .  linesearch - linesearch context
1444f40b411bSPeter Brune 
1445f40b411bSPeter Brune    Output Parameters:
14466232e825SPeter Brune +  X - Solution vector
14476232e825SPeter Brune .  F - Function vector
14486232e825SPeter Brune .  Y - Search direction vector
14496232e825SPeter Brune .  W - Solution work vector
14506232e825SPeter Brune -  G - Function work vector
14516232e825SPeter Brune 
14527bba9028SPeter Brune    Notes:
14537bba9028SPeter Brune    At the beginning of a line search application, X should contain a
14546232e825SPeter Brune    solution and the vector F the function computed at X.  At the end of the
14556232e825SPeter Brune    line search application, X should contain the new solution, and F the
14566232e825SPeter Brune    function evaluated at the new solution.
1457f40b411bSPeter Brune 
14582a7a6963SBarry Smith    These vectors are owned by the SNESLineSearch and should not be destroyed by the caller
14592a7a6963SBarry Smith 
146078bcb3b5SPeter Brune    Level: advanced
1461f40b411bSPeter Brune 
1462f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1463f40b411bSPeter Brune @*/
1464bf388a1fSBarry Smith PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G)
1465bf388a1fSBarry Smith {
14666a388c36SPeter Brune   PetscFunctionBegin;
1467f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
14686a388c36SPeter Brune   if (X) {
14696a388c36SPeter Brune     PetscValidPointer(X, 2);
14706a388c36SPeter Brune     *X = linesearch->vec_sol;
14716a388c36SPeter Brune   }
14726a388c36SPeter Brune   if (F) {
14736a388c36SPeter Brune     PetscValidPointer(F, 3);
14746a388c36SPeter Brune     *F = linesearch->vec_func;
14756a388c36SPeter Brune   }
14766a388c36SPeter Brune   if (Y) {
14776a388c36SPeter Brune     PetscValidPointer(Y, 4);
14786a388c36SPeter Brune     *Y = linesearch->vec_update;
14796a388c36SPeter Brune   }
14806a388c36SPeter Brune   if (W) {
14816a388c36SPeter Brune     PetscValidPointer(W, 5);
14826a388c36SPeter Brune     *W = linesearch->vec_sol_new;
14836a388c36SPeter Brune   }
14846a388c36SPeter Brune   if (G) {
14856a388c36SPeter Brune     PetscValidPointer(G, 6);
14866a388c36SPeter Brune     *G = linesearch->vec_func_new;
14876a388c36SPeter Brune   }
14886a388c36SPeter Brune   PetscFunctionReturn(0);
14896a388c36SPeter Brune }
14906a388c36SPeter Brune 
1491f40b411bSPeter Brune /*@
1492f1c6b773SPeter Brune    SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context
1493f40b411bSPeter Brune 
1494f40b411bSPeter Brune    Input Parameters:
149578bcb3b5SPeter Brune +  linesearch - linesearch context
14966232e825SPeter Brune .  X - Solution vector
14976232e825SPeter Brune .  F - Function vector
14986232e825SPeter Brune .  Y - Search direction vector
14996232e825SPeter Brune .  W - Solution work vector
15006232e825SPeter Brune -  G - Function work vector
1501f40b411bSPeter Brune 
150278bcb3b5SPeter Brune    Level: advanced
1503f40b411bSPeter Brune 
1504f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms(), SNESLineSearchGetVecs()
1505f40b411bSPeter Brune @*/
1506bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G)
1507bf388a1fSBarry Smith {
15086a388c36SPeter Brune   PetscFunctionBegin;
1509f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
15106a388c36SPeter Brune   if (X) {
15116a388c36SPeter Brune     PetscValidHeaderSpecific(X,VEC_CLASSID,2);
15126a388c36SPeter Brune     linesearch->vec_sol = X;
15136a388c36SPeter Brune   }
15146a388c36SPeter Brune   if (F) {
15156a388c36SPeter Brune     PetscValidHeaderSpecific(F,VEC_CLASSID,3);
15166a388c36SPeter Brune     linesearch->vec_func = F;
15176a388c36SPeter Brune   }
15186a388c36SPeter Brune   if (Y) {
15196a388c36SPeter Brune     PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
15206a388c36SPeter Brune     linesearch->vec_update = Y;
15216a388c36SPeter Brune   }
15226a388c36SPeter Brune   if (W) {
15236a388c36SPeter Brune     PetscValidHeaderSpecific(W,VEC_CLASSID,5);
15246a388c36SPeter Brune     linesearch->vec_sol_new = W;
15256a388c36SPeter Brune   }
15266a388c36SPeter Brune   if (G) {
15276a388c36SPeter Brune     PetscValidHeaderSpecific(G,VEC_CLASSID,6);
15286a388c36SPeter Brune     linesearch->vec_func_new = G;
15296a388c36SPeter Brune   }
15306a388c36SPeter Brune   PetscFunctionReturn(0);
15316a388c36SPeter Brune }
15326a388c36SPeter Brune 
1533e7058c64SPeter Brune /*@C
1534f1c6b773SPeter Brune    SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1535e7058c64SPeter Brune    SNES options in the database.
1536e7058c64SPeter Brune 
1537cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
1538e7058c64SPeter Brune 
1539e7058c64SPeter Brune    Input Parameters:
1540e7058c64SPeter Brune +  snes - the SNES context
1541e7058c64SPeter Brune -  prefix - the prefix to prepend to all option names
1542e7058c64SPeter Brune 
1543e7058c64SPeter Brune    Notes:
1544e7058c64SPeter Brune    A hyphen (-) must NOT be given at the beginning of the prefix name.
1545e7058c64SPeter Brune    The first character of all runtime options is AUTOMATICALLY the hyphen.
1546e7058c64SPeter Brune 
1547e7058c64SPeter Brune    Level: advanced
1548e7058c64SPeter Brune 
1549f1c6b773SPeter Brune .keywords: SNESLineSearch, append, options, prefix, database
1550e7058c64SPeter Brune 
1551e7058c64SPeter Brune .seealso: SNESGetOptionsPrefix()
1552e7058c64SPeter Brune @*/
1553f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[])
1554e7058c64SPeter Brune {
1555e7058c64SPeter Brune   PetscErrorCode ierr;
1556e7058c64SPeter Brune 
1557e7058c64SPeter Brune   PetscFunctionBegin;
1558f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1559e7058c64SPeter Brune   ierr = PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1560e7058c64SPeter Brune   PetscFunctionReturn(0);
1561e7058c64SPeter Brune }
1562e7058c64SPeter Brune 
1563e7058c64SPeter Brune /*@C
1564f1c6b773SPeter Brune    SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all
1565f1c6b773SPeter Brune    SNESLineSearch options in the database.
1566e7058c64SPeter Brune 
1567e7058c64SPeter Brune    Not Collective
1568e7058c64SPeter Brune 
1569e7058c64SPeter Brune    Input Parameter:
1570cd7522eaSPeter Brune .  linesearch - the SNESLineSearch context
1571e7058c64SPeter Brune 
1572e7058c64SPeter Brune    Output Parameter:
1573e7058c64SPeter Brune .  prefix - pointer to the prefix string used
1574e7058c64SPeter Brune 
15758e557f58SPeter Brune    Notes:
15768e557f58SPeter Brune    On the fortran side, the user should pass in a string 'prefix' of
1577e7058c64SPeter Brune    sufficient length to hold the prefix.
1578e7058c64SPeter Brune 
1579e7058c64SPeter Brune    Level: advanced
1580e7058c64SPeter Brune 
1581f1c6b773SPeter Brune .keywords: SNESLineSearch, get, options, prefix, database
1582e7058c64SPeter Brune 
1583e7058c64SPeter Brune .seealso: SNESAppendOptionsPrefix()
1584e7058c64SPeter Brune @*/
1585f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[])
1586e7058c64SPeter Brune {
1587e7058c64SPeter Brune   PetscErrorCode ierr;
1588e7058c64SPeter Brune 
1589e7058c64SPeter Brune   PetscFunctionBegin;
1590f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1591e7058c64SPeter Brune   ierr = PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1592e7058c64SPeter Brune   PetscFunctionReturn(0);
1593e7058c64SPeter Brune }
1594bf7f4e0aSPeter Brune 
15958d359177SBarry Smith /*@C
1596fa0ddf94SBarry Smith    SNESLineSearchSetWorkVecs - Gets work vectors for the line search.
1597f40b411bSPeter Brune 
1598f40b411bSPeter Brune    Input Parameter:
1599f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
1600f40b411bSPeter Brune -  nwork - the number of work vectors
1601f40b411bSPeter Brune 
1602f40b411bSPeter Brune    Level: developer
1603f40b411bSPeter Brune 
1604f1c6b773SPeter Brune .keywords: SNESLineSearch, work, vector
1605f40b411bSPeter Brune 
1606fa0ddf94SBarry Smith .seealso: SNESSetWorkVecs()
1607f40b411bSPeter Brune @*/
1608fa0ddf94SBarry Smith PetscErrorCode  SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork)
1609bf7f4e0aSPeter Brune {
1610bf7f4e0aSPeter Brune   PetscErrorCode ierr;
1611bf388a1fSBarry Smith 
1612bf7f4e0aSPeter Brune   PetscFunctionBegin;
1613bf7f4e0aSPeter Brune   if (linesearch->vec_sol) {
1614bf7f4e0aSPeter Brune     ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr);
16158d359177SBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
1616bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1617bf7f4e0aSPeter Brune }
1618bf7f4e0aSPeter Brune 
1619f40b411bSPeter Brune /*@
1620422a814eSBarry Smith    SNESLineSearchGetReason - Gets the success/failure status of the last line search application
1621f40b411bSPeter Brune 
1622f40b411bSPeter Brune    Input Parameters:
162378bcb3b5SPeter Brune .  linesearch - linesearch context
1624f40b411bSPeter Brune 
1625f40b411bSPeter Brune    Output Parameters:
1626422a814eSBarry Smith .  result - The success or failure status
1627f40b411bSPeter Brune 
1628cd7522eaSPeter Brune    Notes:
1629c98378a5SBarry Smith    This is typically called after SNESLineSearchApply() in order to determine if the line-search failed
1630cd7522eaSPeter Brune    (and set the SNES convergence accordingly).
1631cd7522eaSPeter Brune 
1632f40b411bSPeter Brune    Level: intermediate
1633f40b411bSPeter Brune 
1634422a814eSBarry Smith .seealso: SNESLineSearchSetReason(), SNESLineSearchReason
1635f40b411bSPeter Brune @*/
1636422a814eSBarry Smith PetscErrorCode  SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result)
1637bf7f4e0aSPeter Brune {
1638bf7f4e0aSPeter Brune   PetscFunctionBegin;
1639f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1640422a814eSBarry Smith   PetscValidPointer(result, 2);
1641422a814eSBarry Smith   *result = linesearch->result;
1642bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1643bf7f4e0aSPeter Brune }
1644bf7f4e0aSPeter Brune 
1645f40b411bSPeter Brune /*@
1646422a814eSBarry Smith    SNESLineSearchSetReason - Sets the success/failure status of the last line search application
1647f40b411bSPeter Brune 
1648f40b411bSPeter Brune    Input Parameters:
164978bcb3b5SPeter Brune +  linesearch - linesearch context
1650422a814eSBarry Smith -  result - The success or failure status
1651f40b411bSPeter Brune 
1652cd7522eaSPeter Brune    Notes:
1653cd7522eaSPeter Brune    This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set
1654cd7522eaSPeter Brune    the success or failure of the line search method.
1655cd7522eaSPeter Brune 
165678bcb3b5SPeter Brune    Level: developer
1657f40b411bSPeter Brune 
1658422a814eSBarry Smith .seealso: SNESLineSearchGetSResult()
1659f40b411bSPeter Brune @*/
1660422a814eSBarry Smith PetscErrorCode  SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result)
16616a388c36SPeter Brune {
16626a388c36SPeter Brune   PetscFunctionBegin;
16635fd66863SKarl Rupp   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1664422a814eSBarry Smith   linesearch->result = result;
16656a388c36SPeter Brune   PetscFunctionReturn(0);
16666a388c36SPeter Brune }
16676a388c36SPeter Brune 
16689bd66eb0SPeter Brune /*@C
1669f1c6b773SPeter Brune    SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.
16709bd66eb0SPeter Brune 
16719bd66eb0SPeter Brune    Input Parameters:
16729bd66eb0SPeter Brune +  snes - nonlinear context obtained from SNESCreate()
16739bd66eb0SPeter Brune .  projectfunc - function for projecting the function to the bounds
16749bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
16759bd66eb0SPeter Brune 
16769bd66eb0SPeter Brune    Logically Collective on SNES
16779bd66eb0SPeter Brune 
16789bd66eb0SPeter Brune    Calling sequence of projectfunc:
16799bd66eb0SPeter Brune .vb
16809bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X)
16819bd66eb0SPeter Brune .ve
16829bd66eb0SPeter Brune 
16839bd66eb0SPeter Brune     Input parameters for projectfunc:
16849bd66eb0SPeter Brune +   snes - nonlinear context
16859bd66eb0SPeter Brune -   X - current solution
16869bd66eb0SPeter Brune 
1687cd7522eaSPeter Brune     Output parameters for projectfunc:
16889bd66eb0SPeter Brune .   X - Projected solution
16899bd66eb0SPeter Brune 
16909bd66eb0SPeter Brune    Calling sequence of normfunc:
16919bd66eb0SPeter Brune .vb
16929bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm)
16939bd66eb0SPeter Brune .ve
16949bd66eb0SPeter Brune 
1695cd7522eaSPeter Brune     Input parameters for normfunc:
16969bd66eb0SPeter Brune +   snes - nonlinear context
16979bd66eb0SPeter Brune .   X - current solution
16989bd66eb0SPeter Brune -   F - current residual
16999bd66eb0SPeter Brune 
1700cd7522eaSPeter Brune     Output parameters for normfunc:
17019bd66eb0SPeter Brune .   fnorm - VI-specific norm of the function
17029bd66eb0SPeter Brune 
1703cd7522eaSPeter Brune     Notes:
1704cd7522eaSPeter Brune     The VI solvers require projection of the solution to the feasible set.  projectfunc should implement this.
1705cd7522eaSPeter Brune 
1706cd7522eaSPeter Brune     The VI solvers require special evaluation of the function norm such that the norm is only calculated
1707cd7522eaSPeter Brune     on the inactive set.  This should be implemented by normfunc.
17089bd66eb0SPeter Brune 
17099bd66eb0SPeter Brune     Level: developer
17109bd66eb0SPeter Brune 
17119bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, set, line search
17129bd66eb0SPeter Brune 
1713f1c6b773SPeter Brune .seealso: SNESLineSearchGetVIFunctions(), SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck()
17149bd66eb0SPeter Brune @*/
171525acbd8eSLisandro Dalcin PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc)
17169bd66eb0SPeter Brune {
17179bd66eb0SPeter Brune   PetscFunctionBegin;
1718f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
17199bd66eb0SPeter Brune   if (projectfunc) linesearch->ops->viproject = projectfunc;
17209bd66eb0SPeter Brune   if (normfunc) linesearch->ops->vinorm = normfunc;
17219bd66eb0SPeter Brune   PetscFunctionReturn(0);
17229bd66eb0SPeter Brune }
17239bd66eb0SPeter Brune 
17249bd66eb0SPeter Brune /*@C
1725f1c6b773SPeter Brune    SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.
17269bd66eb0SPeter Brune 
17279bd66eb0SPeter Brune    Input Parameters:
1728907376e6SBarry Smith .  linesearch - the line search context, obtain with SNESGetLineSearch()
17299bd66eb0SPeter Brune 
17309bd66eb0SPeter Brune    Output Parameters:
17319bd66eb0SPeter Brune +  projectfunc - function for projecting the function to the bounds
17329bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
17339bd66eb0SPeter Brune 
17349bd66eb0SPeter Brune    Logically Collective on SNES
17359bd66eb0SPeter Brune 
17369bd66eb0SPeter Brune     Level: developer
17379bd66eb0SPeter Brune 
17389bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, get, line search
17399bd66eb0SPeter Brune 
1740f1c6b773SPeter Brune .seealso: SNESLineSearchSetVIFunctions(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
17419bd66eb0SPeter Brune @*/
174225acbd8eSLisandro Dalcin PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc)
17439bd66eb0SPeter Brune {
17449bd66eb0SPeter Brune   PetscFunctionBegin;
17459bd66eb0SPeter Brune   if (projectfunc) *projectfunc = linesearch->ops->viproject;
17469bd66eb0SPeter Brune   if (normfunc) *normfunc = linesearch->ops->vinorm;
17479bd66eb0SPeter Brune   PetscFunctionReturn(0);
17489bd66eb0SPeter Brune }
17499bd66eb0SPeter Brune 
1750bf7f4e0aSPeter Brune /*@C
17511c84c290SBarry Smith   SNESLineSearchRegister - See SNESLineSearchRegister()
1752bf7f4e0aSPeter Brune 
1753bf7f4e0aSPeter Brune   Level: advanced
1754bf7f4e0aSPeter Brune @*/
1755bdf89e91SBarry Smith PetscErrorCode  SNESLineSearchRegister(const char sname[],PetscErrorCode (*function)(SNESLineSearch))
1756bf7f4e0aSPeter Brune {
1757bf7f4e0aSPeter Brune   PetscErrorCode ierr;
1758bf7f4e0aSPeter Brune 
1759bf7f4e0aSPeter Brune   PetscFunctionBegin;
17601d36bdfdSBarry Smith   ierr = SNESInitializePackage();CHKERRQ(ierr);
1761a240a19fSJed Brown   ierr = PetscFunctionListAdd(&SNESLineSearchList,sname,function);CHKERRQ(ierr);
1762bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1763bf7f4e0aSPeter Brune }
1764