xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision 064a246e8b5c1f87897a54b4a9ec05181ea08258)
1af0996ceSBarry Smith #include <petsc/private/linesearchimpl.h> /*I "petscsnes.h" I*/
2bf7f4e0aSPeter Brune 
3f1c6b773SPeter Brune PetscBool         SNESLineSearchRegisterAllCalled = PETSC_FALSE;
40298fd71SBarry Smith PetscFunctionList SNESLineSearchList              = NULL;
5bf7f4e0aSPeter Brune 
6f1c6b773SPeter Brune PetscClassId  SNESLINESEARCH_CLASSID;
757a83faaSBarry Smith PetscLogEvent SNESLINESEARCH_Apply;
8bf7f4e0aSPeter Brune 
9dcf2fd19SBarry Smith /*@
10dcf2fd19SBarry Smith    SNESLineSearchMonitorCancel - Clears all the monitor functions for a SNESLineSearch object.
11dcf2fd19SBarry Smith 
12dcf2fd19SBarry Smith    Logically Collective on SNESLineSearch
13dcf2fd19SBarry Smith 
14dcf2fd19SBarry Smith    Input Parameters:
15dcf2fd19SBarry Smith .  ls - the SNESLineSearch context
16dcf2fd19SBarry Smith 
17dcf2fd19SBarry Smith    Options Database Key:
18dcf2fd19SBarry Smith .  -snes_linesearch_monitor_cancel - cancels all monitors that have been hardwired
19dcf2fd19SBarry Smith     into a code by calls to SNESLineSearchMonitorSet(), but does not cancel those
20dcf2fd19SBarry Smith     set via the options database
21dcf2fd19SBarry Smith 
22dcf2fd19SBarry Smith    Notes:
23dcf2fd19SBarry Smith    There is no way to clear one specific monitor from a SNESLineSearch object.
24dcf2fd19SBarry Smith 
25dcf2fd19SBarry Smith    This does not clear the monitor set with SNESLineSearchSetDefaultMonitor() use SNESLineSearchSetDefaultMonitor(ls,NULL) to cancel
26dcf2fd19SBarry Smith    that one.
27dcf2fd19SBarry Smith 
28dcf2fd19SBarry Smith    Level: intermediate
29dcf2fd19SBarry Smith 
3084238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchMonitorDefault(), SNESLineSearchMonitorSet()
31dcf2fd19SBarry Smith @*/
32dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchMonitorCancel(SNESLineSearch ls)
33dcf2fd19SBarry Smith {
34dcf2fd19SBarry Smith   PetscErrorCode ierr;
35dcf2fd19SBarry Smith   PetscInt       i;
36dcf2fd19SBarry Smith 
37dcf2fd19SBarry Smith   PetscFunctionBegin;
38dcf2fd19SBarry Smith   PetscValidHeaderSpecific(ls,SNESLINESEARCH_CLASSID,1);
39dcf2fd19SBarry Smith   for (i=0; i<ls->numbermonitors; i++) {
40dcf2fd19SBarry Smith     if (ls->monitordestroy[i]) {
41dcf2fd19SBarry Smith       ierr = (*ls->monitordestroy[i])(&ls->monitorcontext[i]);CHKERRQ(ierr);
42dcf2fd19SBarry Smith     }
43dcf2fd19SBarry Smith   }
44dcf2fd19SBarry Smith   ls->numbermonitors = 0;
45dcf2fd19SBarry Smith   PetscFunctionReturn(0);
46dcf2fd19SBarry Smith }
47dcf2fd19SBarry Smith 
48dcf2fd19SBarry Smith /*@
49dcf2fd19SBarry Smith    SNESLineSearchMonitor - runs the user provided monitor routines, if they exist
50dcf2fd19SBarry Smith 
51dcf2fd19SBarry Smith    Collective on SNES
52dcf2fd19SBarry Smith 
53dcf2fd19SBarry Smith    Input Parameters:
54dcf2fd19SBarry Smith .  ls - the linesearch object
55dcf2fd19SBarry Smith 
56dcf2fd19SBarry Smith    Notes:
57dcf2fd19SBarry Smith    This routine is called by the SNES implementations.
58dcf2fd19SBarry Smith    It does not typically need to be called by the user.
59dcf2fd19SBarry Smith 
60dcf2fd19SBarry Smith    Level: developer
61dcf2fd19SBarry Smith 
6284238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchMonitorSet()
63dcf2fd19SBarry Smith @*/
64dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchMonitor(SNESLineSearch ls)
65dcf2fd19SBarry Smith {
66dcf2fd19SBarry Smith   PetscErrorCode ierr;
67dcf2fd19SBarry Smith   PetscInt       i,n = ls->numbermonitors;
68dcf2fd19SBarry Smith 
69dcf2fd19SBarry Smith   PetscFunctionBegin;
70dcf2fd19SBarry Smith   for (i=0; i<n; i++) {
71dcf2fd19SBarry Smith     ierr = (*ls->monitorftns[i])(ls,ls->monitorcontext[i]);CHKERRQ(ierr);
72dcf2fd19SBarry Smith   }
73dcf2fd19SBarry Smith   PetscFunctionReturn(0);
74dcf2fd19SBarry Smith }
75dcf2fd19SBarry Smith 
76dcf2fd19SBarry Smith /*@C
77dcf2fd19SBarry Smith    SNESLineSearchMonitorSet - Sets an ADDITIONAL function that is to be used at every
78dcf2fd19SBarry Smith    iteration of the nonlinear solver to display the iteration's
79dcf2fd19SBarry Smith    progress.
80dcf2fd19SBarry Smith 
81dcf2fd19SBarry Smith    Logically Collective on SNESLineSearch
82dcf2fd19SBarry Smith 
83dcf2fd19SBarry Smith    Input Parameters:
84dcf2fd19SBarry Smith +  ls - the SNESLineSearch context
85dcf2fd19SBarry Smith .  f - the monitor function
86dcf2fd19SBarry Smith .  mctx - [optional] user-defined context for private data for the
87dcf2fd19SBarry Smith           monitor routine (use NULL if no context is desired)
88dcf2fd19SBarry Smith -  monitordestroy - [optional] routine that frees monitor context
89dcf2fd19SBarry Smith           (may be NULL)
90dcf2fd19SBarry Smith 
91dcf2fd19SBarry Smith    Notes:
92dcf2fd19SBarry Smith    Several different monitoring routines may be set by calling
93dcf2fd19SBarry Smith    SNESLineSearchMonitorSet() multiple times; all will be called in the
94dcf2fd19SBarry Smith    order in which they were set.
95dcf2fd19SBarry Smith 
9695452b02SPatrick Sanan    Fortran Notes:
9795452b02SPatrick Sanan     Only a single monitor function can be set for each SNESLineSearch object
98dcf2fd19SBarry Smith 
99dcf2fd19SBarry Smith    Level: intermediate
100dcf2fd19SBarry Smith 
10184238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchMonitorDefault(), SNESLineSearchMonitorCancel()
102dcf2fd19SBarry Smith @*/
103dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchMonitorSet(SNESLineSearch ls,PetscErrorCode (*f)(SNESLineSearch,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
104dcf2fd19SBarry Smith {
10578064530SBarry Smith   PetscErrorCode ierr;
10678064530SBarry Smith   PetscInt       i;
10778064530SBarry Smith   PetscBool      identical;
10878064530SBarry Smith 
109dcf2fd19SBarry Smith   PetscFunctionBegin;
110dcf2fd19SBarry Smith   PetscValidHeaderSpecific(ls,SNESLINESEARCH_CLASSID,1);
11178064530SBarry Smith   for (i=0; i<ls->numbermonitors;i++) {
11278064530SBarry Smith     ierr = PetscMonitorCompare((PetscErrorCode (*)(void))f,mctx,monitordestroy,(PetscErrorCode (*)(void))ls->monitorftns[i],ls->monitorcontext[i],ls->monitordestroy[i],&identical);CHKERRQ(ierr);
11378064530SBarry Smith     if (identical) PetscFunctionReturn(0);
11478064530SBarry Smith   }
115dcf2fd19SBarry Smith   if (ls->numbermonitors >= MAXSNESLSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
116dcf2fd19SBarry Smith   ls->monitorftns[ls->numbermonitors]          = f;
117dcf2fd19SBarry Smith   ls->monitordestroy[ls->numbermonitors]   = monitordestroy;
118dcf2fd19SBarry Smith   ls->monitorcontext[ls->numbermonitors++] = (void*)mctx;
119dcf2fd19SBarry Smith   PetscFunctionReturn(0);
120dcf2fd19SBarry Smith }
121dcf2fd19SBarry Smith 
122dcf2fd19SBarry Smith /*@C
123dcf2fd19SBarry Smith    SNESLineSearchMonitorSolutionUpdate - Monitors each update a new function value the linesearch tries
124dcf2fd19SBarry Smith 
125dcf2fd19SBarry Smith    Collective on SNESLineSearch
126dcf2fd19SBarry Smith 
127dcf2fd19SBarry Smith    Input Parameters:
128dcf2fd19SBarry Smith +  ls - the SNES linesearch object
129d12e167eSBarry Smith -  vf - the context for the monitor, in this case it is an ASCII PetscViewer and format
130dcf2fd19SBarry Smith 
131dcf2fd19SBarry Smith    Level: intermediate
132dcf2fd19SBarry Smith 
13384238204SBarry Smith .seealso: SNESGetLineSearch(), SNESMonitorSet(), SNESMonitorSolution()
134dcf2fd19SBarry Smith @*/
135d12e167eSBarry Smith PetscErrorCode  SNESLineSearchMonitorSolutionUpdate(SNESLineSearch ls,PetscViewerAndFormat *vf)
136dcf2fd19SBarry Smith {
137dcf2fd19SBarry Smith   PetscErrorCode ierr;
138d12e167eSBarry Smith   PetscViewer    viewer = vf->viewer;
139dcf2fd19SBarry Smith   Vec            Y,W,G;
140dcf2fd19SBarry Smith 
141dcf2fd19SBarry Smith   PetscFunctionBegin;
142dcf2fd19SBarry Smith   ierr = SNESLineSearchGetVecs(ls,NULL,NULL,&Y,&W,&G);CHKERRQ(ierr);
143d12e167eSBarry Smith   ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr);
144dcf2fd19SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"LineSearch attempted update to solution \n");CHKERRQ(ierr);
145dcf2fd19SBarry Smith   ierr = VecView(Y,viewer);CHKERRQ(ierr);
146dcf2fd19SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"LineSearch attempted new solution \n");CHKERRQ(ierr);
147dcf2fd19SBarry Smith   ierr = VecView(W,viewer);CHKERRQ(ierr);
148dcf2fd19SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"LineSearch attempted updated function value\n");CHKERRQ(ierr);
149dcf2fd19SBarry Smith   ierr = VecView(G,viewer);CHKERRQ(ierr);
150d12e167eSBarry Smith   ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
151dcf2fd19SBarry Smith   PetscFunctionReturn(0);
152dcf2fd19SBarry Smith }
153dcf2fd19SBarry Smith 
154f40b411bSPeter Brune /*@
155cd7522eaSPeter Brune    SNESLineSearchCreate - Creates the line search context.
156f40b411bSPeter Brune 
157cd7522eaSPeter Brune    Logically Collective on Comm
158f40b411bSPeter Brune 
159f40b411bSPeter Brune    Input Parameters:
160cd7522eaSPeter Brune .  comm - MPI communicator for the line search (typically from the associated SNES context).
161f40b411bSPeter Brune 
162f40b411bSPeter Brune    Output Parameters:
1638e557f58SPeter Brune .  outlinesearch - the new linesearch context
164f40b411bSPeter Brune 
165162e0bf5SPeter Brune    Level: developer
166162e0bf5SPeter Brune 
167162e0bf5SPeter Brune    Notes:
168162e0bf5SPeter Brune    The preferred calling sequence for users is to use SNESGetLineSearch() to acquire the SNESLineSearch instance
169162e0bf5SPeter Brune    already associated with the SNES.  This function is for developer use.
170f40b411bSPeter Brune 
171162e0bf5SPeter Brune .seealso: LineSearchDestroy(), SNESGetLineSearch()
172f40b411bSPeter Brune @*/
173f40b411bSPeter Brune 
174bf388a1fSBarry Smith PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch)
175bf388a1fSBarry Smith {
176bf7f4e0aSPeter Brune   PetscErrorCode ierr;
177f1c6b773SPeter Brune   SNESLineSearch linesearch;
178bf388a1fSBarry Smith 
179bf7f4e0aSPeter Brune   PetscFunctionBegin;
180ea5d4fccSPeter Brune   PetscValidPointer(outlinesearch,2);
181f34a81feSMatthew G. Knepley   ierr = SNESInitializePackage();CHKERRQ(ierr);
1820298fd71SBarry Smith   *outlinesearch = NULL;
183f5af7f23SKarl Rupp 
18473107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(linesearch,SNESLINESEARCH_CLASSID, "SNESLineSearch","Linesearch","SNESLineSearch",comm,SNESLineSearchDestroy,SNESLineSearchView);CHKERRQ(ierr);
185bf7f4e0aSPeter Brune 
1860298fd71SBarry Smith   linesearch->vec_sol_new  = NULL;
1870298fd71SBarry Smith   linesearch->vec_func_new = NULL;
1880298fd71SBarry Smith   linesearch->vec_sol      = NULL;
1890298fd71SBarry Smith   linesearch->vec_func     = NULL;
1900298fd71SBarry Smith   linesearch->vec_update   = NULL;
1919bd66eb0SPeter Brune 
192bf7f4e0aSPeter Brune   linesearch->lambda       = 1.0;
193bf7f4e0aSPeter Brune   linesearch->fnorm        = 1.0;
194bf7f4e0aSPeter Brune   linesearch->ynorm        = 1.0;
195bf7f4e0aSPeter Brune   linesearch->xnorm        = 1.0;
196422a814eSBarry Smith   linesearch->result       = SNES_LINESEARCH_SUCCEEDED;
197bf7f4e0aSPeter Brune   linesearch->norms        = PETSC_TRUE;
198bf7f4e0aSPeter Brune   linesearch->keeplambda   = PETSC_FALSE;
199bf7f4e0aSPeter Brune   linesearch->damping      = 1.0;
200bf7f4e0aSPeter Brune   linesearch->maxstep      = 1e8;
201bf7f4e0aSPeter Brune   linesearch->steptol      = 1e-12;
202516fe3c3SPeter Brune   linesearch->rtol         = 1e-8;
203516fe3c3SPeter Brune   linesearch->atol         = 1e-15;
204516fe3c3SPeter Brune   linesearch->ltol         = 1e-8;
2050298fd71SBarry Smith   linesearch->precheckctx  = NULL;
2060298fd71SBarry Smith   linesearch->postcheckctx = NULL;
20759405d5eSPeter Brune   linesearch->max_its      = 1;
208bf7f4e0aSPeter Brune   linesearch->setupcalled  = PETSC_FALSE;
2093add74b1SBarry Smith   linesearch->monitor      = NULL;
210bf7f4e0aSPeter Brune   *outlinesearch           = linesearch;
211bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
212bf7f4e0aSPeter Brune }
213bf7f4e0aSPeter Brune 
214f40b411bSPeter Brune /*@
21578bcb3b5SPeter Brune    SNESLineSearchSetUp - Prepares the line search for being applied by allocating
21678bcb3b5SPeter Brune    any required vectors.
217f40b411bSPeter Brune 
218cd7522eaSPeter Brune    Collective on SNESLineSearch
219f40b411bSPeter Brune 
220f40b411bSPeter Brune    Input Parameters:
221f40b411bSPeter Brune .  linesearch - The LineSearch instance.
222f40b411bSPeter Brune 
223cd7522eaSPeter Brune    Notes:
224f190f2fcSBarry Smith    For most cases, this needn't be called by users or outside of SNESLineSearchApply().
225cd7522eaSPeter Brune    The only current case where this is called outside of this is for the VI
22678bcb3b5SPeter Brune    solvers, which modify the solution and work vectors before the first call
227cd7522eaSPeter Brune    of SNESLineSearchApply, requiring the SNESLineSearch work vectors to be
228cd7522eaSPeter Brune    allocated upfront.
229cd7522eaSPeter Brune 
23078bcb3b5SPeter Brune    Level: advanced
231f40b411bSPeter Brune 
23284238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchReset()
233f40b411bSPeter Brune @*/
234f40b411bSPeter Brune 
235bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch)
236bf388a1fSBarry Smith {
237bf7f4e0aSPeter Brune   PetscErrorCode ierr;
238bf388a1fSBarry Smith 
239bf7f4e0aSPeter Brune   PetscFunctionBegin;
240bf7f4e0aSPeter Brune   if (!((PetscObject)linesearch)->type_name) {
2411a4f838cSPeter Brune     ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr);
242bf7f4e0aSPeter Brune   }
243bf7f4e0aSPeter Brune   if (!linesearch->setupcalled) {
2449bd66eb0SPeter Brune     if (!linesearch->vec_sol_new) {
245bf7f4e0aSPeter Brune       ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new);CHKERRQ(ierr);
2469bd66eb0SPeter Brune     }
2479bd66eb0SPeter Brune     if (!linesearch->vec_func_new) {
2489bd66eb0SPeter Brune       ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new);CHKERRQ(ierr);
2499bd66eb0SPeter Brune     }
250bf7f4e0aSPeter Brune     if (linesearch->ops->setup) {
251bf7f4e0aSPeter Brune       ierr = (*linesearch->ops->setup)(linesearch);CHKERRQ(ierr);
252bf7f4e0aSPeter Brune     }
253ed07d7d7SPeter Brune     if (!linesearch->ops->snesfunc) {ierr = SNESLineSearchSetFunction(linesearch,SNESComputeFunction);CHKERRQ(ierr);}
254bf7f4e0aSPeter Brune     linesearch->lambda      = linesearch->damping;
255bf7f4e0aSPeter Brune     linesearch->setupcalled = PETSC_TRUE;
256bf7f4e0aSPeter Brune   }
257bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
258bf7f4e0aSPeter Brune }
259bf7f4e0aSPeter Brune 
260f40b411bSPeter Brune 
261f40b411bSPeter Brune /*@
262f190f2fcSBarry Smith    SNESLineSearchReset - Undoes the SNESLineSearchSetUp() and deletes any Vecs or Mats allocated by the line search.
263f40b411bSPeter Brune 
264f1c6b773SPeter Brune    Collective on SNESLineSearch
265f40b411bSPeter Brune 
266f40b411bSPeter Brune    Input Parameters:
267f40b411bSPeter Brune .  linesearch - The LineSearch instance.
268f40b411bSPeter Brune 
26995452b02SPatrick Sanan    Notes:
27095452b02SPatrick Sanan     Usually only called by SNESReset()
271f190f2fcSBarry Smith 
272f190f2fcSBarry Smith    Level: developer
273f40b411bSPeter Brune 
27484238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchSetUp()
275f40b411bSPeter Brune @*/
276f40b411bSPeter Brune 
277bf388a1fSBarry Smith PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch)
278bf388a1fSBarry Smith {
279bf7f4e0aSPeter Brune   PetscErrorCode ierr;
280bf388a1fSBarry Smith 
281bf7f4e0aSPeter Brune   PetscFunctionBegin;
282f5af7f23SKarl Rupp   if (linesearch->ops->reset) (*linesearch->ops->reset)(linesearch);
283f5af7f23SKarl Rupp 
284bf7f4e0aSPeter Brune   ierr = VecDestroy(&linesearch->vec_sol_new);CHKERRQ(ierr);
285bf7f4e0aSPeter Brune   ierr = VecDestroy(&linesearch->vec_func_new);CHKERRQ(ierr);
286bf7f4e0aSPeter Brune 
287bf7f4e0aSPeter Brune   ierr = VecDestroyVecs(linesearch->nwork, &linesearch->work);CHKERRQ(ierr);
288f5af7f23SKarl Rupp 
289bf7f4e0aSPeter Brune   linesearch->nwork       = 0;
290bf7f4e0aSPeter Brune   linesearch->setupcalled = PETSC_FALSE;
291bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
292bf7f4e0aSPeter Brune }
293bf7f4e0aSPeter Brune 
294ed07d7d7SPeter Brune /*@C
295f190f2fcSBarry Smith    SNESLineSearchSetFunction - Sets the function evaluation used by the SNES line search
296ed07d7d7SPeter Brune 
297ed07d7d7SPeter Brune    Input Parameters:
298ed07d7d7SPeter Brune .  linesearch - the SNESLineSearch context
299f190f2fcSBarry Smith +  func       - function evaluation routine
300ed07d7d7SPeter Brune 
301ed07d7d7SPeter Brune    Level: developer
302ed07d7d7SPeter Brune 
30395452b02SPatrick Sanan    Notes:
30495452b02SPatrick Sanan     This is used internally by PETSc and not called by users
305f190f2fcSBarry Smith 
30684238204SBarry Smith .seealso: SNESGetLineSearch(), SNESSetFunction()
307ed07d7d7SPeter Brune @*/
308ed07d7d7SPeter Brune PetscErrorCode  SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES,Vec,Vec))
309ed07d7d7SPeter Brune {
310ed07d7d7SPeter Brune   PetscFunctionBegin;
311ed07d7d7SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
312ed07d7d7SPeter Brune   linesearch->ops->snesfunc = func;
313ed07d7d7SPeter Brune   PetscFunctionReturn(0);
314ed07d7d7SPeter Brune }
315ed07d7d7SPeter Brune 
31686d74e61SPeter Brune /*@C
317f190f2fcSBarry Smith    SNESLineSearchSetPreCheck - Sets a user function that is called after the initial search direction has been computed but
318df3898eeSBarry Smith          before the line search routine has been applied. Allows the user to adjust the result of (usually a linear solve) that
319f190f2fcSBarry Smith          determined the search direction.
32086d74e61SPeter Brune 
321f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
32286d74e61SPeter Brune 
32386d74e61SPeter Brune    Input Parameters:
324f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
32584238204SBarry Smith .  func - [optional] function evaluation routine, see SNESLineSearchPreCheck() for the calling sequence
326f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
32786d74e61SPeter Brune 
32886d74e61SPeter Brune    Level: intermediate
32986d74e61SPeter Brune 
33084238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
33186d74e61SPeter Brune @*/
332f190f2fcSBarry Smith PetscErrorCode  SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void *ctx)
33386d74e61SPeter Brune {
3349bd66eb0SPeter Brune   PetscFunctionBegin;
335f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
336f190f2fcSBarry Smith   if (func) linesearch->ops->precheck = func;
33786d74e61SPeter Brune   if (ctx) linesearch->precheckctx = ctx;
33886d74e61SPeter Brune   PetscFunctionReturn(0);
33986d74e61SPeter Brune }
34086d74e61SPeter Brune 
34186d74e61SPeter Brune /*@C
34253deb899SBarry Smith    SNESLineSearchGetPreCheck - Gets the pre-check function for the line search routine.
34386d74e61SPeter Brune 
34486d74e61SPeter Brune    Input Parameters:
345f1c6b773SPeter Brune .  linesearch - the SNESLineSearch context
34686d74e61SPeter Brune 
34786d74e61SPeter Brune    Output Parameters:
34884238204SBarry Smith +  func       - [optional] function evaluation routine, see SNESLineSearchPreCheck() for calling sequence
349f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
35086d74e61SPeter Brune 
35186d74e61SPeter Brune    Level: intermediate
35286d74e61SPeter Brune 
35384238204SBarry Smith .seealso: SNESGetLineSearch(), SNESGetLineSearch(), SNESLineSearchPreCheck(), SNESLineSearchGetPostCheck(), SNESLineSearchSetPreCheck(), SNESLineSearchSetPostCheck()
35486d74e61SPeter Brune @*/
355f190f2fcSBarry Smith PetscErrorCode  SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void **ctx)
35686d74e61SPeter Brune {
3579bd66eb0SPeter Brune   PetscFunctionBegin;
358f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
359f190f2fcSBarry Smith   if (func) *func = linesearch->ops->precheck;
36086d74e61SPeter Brune   if (ctx) *ctx = linesearch->precheckctx;
36186d74e61SPeter Brune   PetscFunctionReturn(0);
36286d74e61SPeter Brune }
36386d74e61SPeter Brune 
36486d74e61SPeter Brune /*@C
365f190f2fcSBarry Smith    SNESLineSearchSetPostCheck - Sets a user function that is called after the line search has been applied to determine the step
366f190f2fcSBarry Smith        direction and length. Allows the user a chance to change or override the decision of the line search routine
36786d74e61SPeter Brune 
368f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
36986d74e61SPeter Brune 
37086d74e61SPeter Brune    Input Parameters:
371f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
37284238204SBarry Smith .  func - [optional] function evaluation routine, see SNESLineSearchPostCheck()  for the calling sequence
373f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
37486d74e61SPeter Brune 
37586d74e61SPeter Brune    Level: intermediate
37686d74e61SPeter Brune 
37784238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchPostCheck(), SNESLineSearchSetPreCheck(), SNESLineSearchGetPreCheck(), SNESLineSearchGetPostCheck()
37886d74e61SPeter Brune @*/
379f190f2fcSBarry Smith PetscErrorCode  SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void *ctx)
38086d74e61SPeter Brune {
38186d74e61SPeter Brune   PetscFunctionBegin;
382f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
383f190f2fcSBarry Smith   if (func) linesearch->ops->postcheck = func;
38486d74e61SPeter Brune   if (ctx) linesearch->postcheckctx = ctx;
38586d74e61SPeter Brune   PetscFunctionReturn(0);
38686d74e61SPeter Brune }
38786d74e61SPeter Brune 
38886d74e61SPeter Brune /*@C
389f1c6b773SPeter Brune    SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine.
39086d74e61SPeter Brune 
39186d74e61SPeter Brune    Input Parameters:
392f1c6b773SPeter Brune .  linesearch - the SNESLineSearch context
39386d74e61SPeter Brune 
39486d74e61SPeter Brune    Output Parameters:
39584238204SBarry Smith +  func - [optional] function evaluation routine, see for the calling sequence SNESLineSearchPostCheck()
396f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
39786d74e61SPeter Brune 
39886d74e61SPeter Brune    Level: intermediate
39986d74e61SPeter Brune 
40084238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchPostCheck(), SNESLineSearchSetPreCheck()
40186d74e61SPeter Brune @*/
402f190f2fcSBarry Smith PetscErrorCode  SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void **ctx)
40386d74e61SPeter Brune {
4049bd66eb0SPeter Brune   PetscFunctionBegin;
405f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
406f190f2fcSBarry Smith   if (func) *func = linesearch->ops->postcheck;
40786d74e61SPeter Brune   if (ctx) *ctx = linesearch->postcheckctx;
40886d74e61SPeter Brune   PetscFunctionReturn(0);
40986d74e61SPeter Brune }
41086d74e61SPeter Brune 
411f40b411bSPeter Brune /*@
412f1c6b773SPeter Brune    SNESLineSearchPreCheck - Prepares the line search for being applied.
413f40b411bSPeter Brune 
414cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
415f40b411bSPeter Brune 
416f40b411bSPeter Brune    Input Parameters:
4177b1df9c1SPeter Brune +  linesearch - The linesearch instance.
4187b1df9c1SPeter Brune .  X - The current solution
4197b1df9c1SPeter Brune -  Y - The step direction
420f40b411bSPeter Brune 
421f40b411bSPeter Brune    Output Parameters:
4228e557f58SPeter Brune .  changed - Indicator that the precheck routine has changed anything
423f40b411bSPeter Brune 
424f190f2fcSBarry Smith    Level: developer
425f40b411bSPeter Brune 
42684238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchPostCheck(), SNESLineSearchSetPreCheck(), SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck()
427f40b411bSPeter Brune @*/
4287b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed)
429bf7f4e0aSPeter Brune {
430bf7f4e0aSPeter Brune   PetscErrorCode ierr;
4315fd66863SKarl Rupp 
432bf7f4e0aSPeter Brune   PetscFunctionBegin;
433bf7f4e0aSPeter Brune   *changed = PETSC_FALSE;
4346b2b7091SBarry Smith   if (linesearch->ops->precheck) {
4356b2b7091SBarry Smith     ierr = (*linesearch->ops->precheck)(linesearch, X, Y, changed, linesearch->precheckctx);CHKERRQ(ierr);
43638bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed,4);
437bf7f4e0aSPeter Brune   }
438bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
439bf7f4e0aSPeter Brune }
440bf7f4e0aSPeter Brune 
441f40b411bSPeter Brune /*@
442f1c6b773SPeter Brune    SNESLineSearchPostCheck - Prepares the line search for being applied.
443f40b411bSPeter Brune 
444cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
445f40b411bSPeter Brune 
446f40b411bSPeter Brune    Input Parameters:
4477b1df9c1SPeter Brune +  linesearch - The linesearch context
4487b1df9c1SPeter Brune .  X - The last solution
4497b1df9c1SPeter Brune .  Y - The step direction
4507b1df9c1SPeter Brune -  W - The updated solution, W = X + lambda*Y for some lambda
451f40b411bSPeter Brune 
452f40b411bSPeter Brune    Output Parameters:
45378bcb3b5SPeter Brune +  changed_Y - Indicator if the direction Y has been changed.
45478bcb3b5SPeter Brune -  changed_W - Indicator if the new candidate solution W has been changed.
455f40b411bSPeter Brune 
456f190f2fcSBarry Smith    Level: developer
457f40b411bSPeter Brune 
45884238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck(), SNESLineSearchSetPrecheck(), SNESLineSearchGetPrecheck()
459f40b411bSPeter Brune @*/
4607b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch,Vec X,Vec Y,Vec W,PetscBool *changed_Y,PetscBool *changed_W)
461bf7f4e0aSPeter Brune {
462bf7f4e0aSPeter Brune   PetscErrorCode ierr;
463bf388a1fSBarry Smith 
464bf7f4e0aSPeter Brune   PetscFunctionBegin;
465bf7f4e0aSPeter Brune   *changed_Y = PETSC_FALSE;
466bf7f4e0aSPeter Brune   *changed_W = PETSC_FALSE;
4676b2b7091SBarry Smith   if (linesearch->ops->postcheck) {
4686b2b7091SBarry Smith     ierr = (*linesearch->ops->postcheck)(linesearch,X,Y,W,changed_Y,changed_W,linesearch->postcheckctx);CHKERRQ(ierr);
46938bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed_Y,5);
47038bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed_W,6);
47186d74e61SPeter Brune   }
47286d74e61SPeter Brune   PetscFunctionReturn(0);
47386d74e61SPeter Brune }
47486d74e61SPeter Brune 
47586d74e61SPeter Brune /*@C
47686d74e61SPeter Brune    SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration
47786d74e61SPeter Brune 
478cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
47986d74e61SPeter Brune 
48086d74e61SPeter Brune    Input Arguments:
48186d74e61SPeter Brune +  linesearch - linesearch context
48286d74e61SPeter Brune .  X - base state for this step
48386d74e61SPeter Brune .  Y - initial correction
484907376e6SBarry Smith -  ctx - context for this function
48586d74e61SPeter Brune 
48686d74e61SPeter Brune    Output Arguments:
48786d74e61SPeter Brune +  Y - correction, possibly modified
48886d74e61SPeter Brune -  changed - flag indicating that Y was modified
48986d74e61SPeter Brune 
49086d74e61SPeter Brune    Options Database Key:
491cd7522eaSPeter Brune +  -snes_linesearch_precheck_picard - activate this routine
492cd7522eaSPeter Brune -  -snes_linesearch_precheck_picard_angle - angle
49386d74e61SPeter Brune 
49486d74e61SPeter Brune    Level: advanced
49586d74e61SPeter Brune 
49686d74e61SPeter Brune    Notes:
49786d74e61SPeter Brune    This function should be passed to SNESLineSearchSetPreCheck()
49886d74e61SPeter Brune 
49986d74e61SPeter Brune    The justification for this method involves the linear convergence of a Picard iteration
50086d74e61SPeter Brune    so the Picard linearization should be provided in place of the "Jacobian". This correction
50186d74e61SPeter Brune    is generally not useful when using a Newton linearization.
50286d74e61SPeter Brune 
50386d74e61SPeter Brune    Reference:
50486d74e61SPeter Brune    Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology.
50586d74e61SPeter Brune 
50684238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchSetPreCheck()
50786d74e61SPeter Brune @*/
508f1c6b773SPeter Brune PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed,void *ctx)
50986d74e61SPeter Brune {
51086d74e61SPeter Brune   PetscErrorCode ierr;
51186d74e61SPeter Brune   PetscReal      angle = *(PetscReal*)linesearch->precheckctx;
51286d74e61SPeter Brune   Vec            Ylast;
51386d74e61SPeter Brune   PetscScalar    dot;
51486d74e61SPeter Brune   PetscInt       iter;
51586d74e61SPeter Brune   PetscReal      ynorm,ylastnorm,theta,angle_radians;
51686d74e61SPeter Brune   SNES           snes;
51786d74e61SPeter Brune 
51886d74e61SPeter Brune   PetscFunctionBegin;
519f1c6b773SPeter Brune   ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
52086d74e61SPeter Brune   ierr = PetscObjectQuery((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject*)&Ylast);CHKERRQ(ierr);
52186d74e61SPeter Brune   if (!Ylast) {
52286d74e61SPeter Brune     ierr = VecDuplicate(Y,&Ylast);CHKERRQ(ierr);
52386d74e61SPeter Brune     ierr = PetscObjectCompose((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject)Ylast);CHKERRQ(ierr);
52486d74e61SPeter Brune     ierr = PetscObjectDereference((PetscObject)Ylast);CHKERRQ(ierr);
52586d74e61SPeter Brune   }
52686d74e61SPeter Brune   ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr);
52786d74e61SPeter Brune   if (iter < 2) {
52886d74e61SPeter Brune     ierr     = VecCopy(Y,Ylast);CHKERRQ(ierr);
52986d74e61SPeter Brune     *changed = PETSC_FALSE;
53086d74e61SPeter Brune     PetscFunctionReturn(0);
53186d74e61SPeter Brune   }
53286d74e61SPeter Brune 
53386d74e61SPeter Brune   ierr = VecDot(Y,Ylast,&dot);CHKERRQ(ierr);
53486d74e61SPeter Brune   ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr);
53586d74e61SPeter Brune   ierr = VecNorm(Ylast,NORM_2,&ylastnorm);CHKERRQ(ierr);
53686d74e61SPeter Brune   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
537255453a1SBarry Smith   theta         = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0));
53886d74e61SPeter Brune   angle_radians = angle * PETSC_PI / 180.;
53986d74e61SPeter Brune   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
54086d74e61SPeter Brune     /* Modify the step Y */
54186d74e61SPeter Brune     PetscReal alpha,ydiffnorm;
54286d74e61SPeter Brune     ierr  = VecAXPY(Ylast,-1.0,Y);CHKERRQ(ierr);
54386d74e61SPeter Brune     ierr  = VecNorm(Ylast,NORM_2,&ydiffnorm);CHKERRQ(ierr);
544f85e2ce2SBarry Smith     alpha = (ydiffnorm > .001*ylastnorm) ? ylastnorm / ydiffnorm : 1000.0;
54586d74e61SPeter Brune     ierr  = VecCopy(Y,Ylast);CHKERRQ(ierr);
54686d74e61SPeter Brune     ierr  = VecScale(Y,alpha);CHKERRQ(ierr);
547c69d1a72SBarry 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);
548a47ec85fSBarry Smith     *changed = PETSC_TRUE;
54986d74e61SPeter Brune   } else {
550c69d1a72SBarry 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);
55186d74e61SPeter Brune     ierr     = VecCopy(Y,Ylast);CHKERRQ(ierr);
55286d74e61SPeter Brune     *changed = PETSC_FALSE;
553bf7f4e0aSPeter Brune   }
554bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
555bf7f4e0aSPeter Brune }
556bf7f4e0aSPeter Brune 
557f40b411bSPeter Brune /*@
558cd7522eaSPeter Brune    SNESLineSearchApply - Computes the line-search update.
559f40b411bSPeter Brune 
560f1c6b773SPeter Brune    Collective on SNESLineSearch
561f40b411bSPeter Brune 
562f40b411bSPeter Brune    Input Parameters:
5638e557f58SPeter Brune +  linesearch - The linesearch context
5648e557f58SPeter Brune .  X - The current solution
5658e557f58SPeter Brune .  F - The current function
5668e557f58SPeter Brune .  fnorm - The current norm
5678e557f58SPeter Brune -  Y - The search direction
568f40b411bSPeter Brune 
569f40b411bSPeter Brune    Output Parameters:
5708e557f58SPeter Brune +  X - The new solution
5718e557f58SPeter Brune .  F - The new function
5728e557f58SPeter Brune -  fnorm - The new function norm
573f40b411bSPeter Brune 
574cd7522eaSPeter Brune    Options Database Keys:
575d4c6564cSPatrick Farrell + -snes_linesearch_type - basic, bt, l2, cp, nleqerr, shell
576dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches
5771fe24845SBarry Smith . -snes_linesearch_damping - The linesearch damping parameter, default is 1.0 (no damping)
5781fe24845SBarry Smith . -snes_linesearch_norms   - Turn on/off the linesearch norms computation (SNESLineSearchSetComputeNorms())
5793c7d6663SPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess
5803c7d6663SPeter Brune - -snes_linesearch_max_it - The number of iterations for iterative line searches
581cd7522eaSPeter Brune 
582cd7522eaSPeter Brune    Notes:
583cd7522eaSPeter Brune    This is typically called from within a SNESSolve() implementation in order to
584cd7522eaSPeter Brune    help with convergence of the nonlinear method.  Various SNES types use line searches
585cd7522eaSPeter Brune    in different ways, but the overarching theme is that a line search is used to determine
586cd7522eaSPeter Brune    an optimal damping parameter of a step at each iteration of the method.  Each
58784238204SBarry Smith    application of the line search may invoke SNESComputeFunction() several times, and
588cd7522eaSPeter Brune    therefore may be fairly expensive.
589cd7522eaSPeter Brune 
590f40b411bSPeter Brune    Level: Intermediate
591f40b411bSPeter Brune 
59284238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchCreate(), SNESLineSearchPreCheck(), SNESLineSearchPostCheck(), SNESSolve(), SNESComputeFunction(), SNESLineSearchSetComputeNorms(),
5931fe24845SBarry Smith           SNESLineSearchType, SNESLineSearchSetType()
594f40b411bSPeter Brune @*/
595bf388a1fSBarry Smith PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y)
596bf388a1fSBarry Smith {
597bf7f4e0aSPeter Brune   PetscErrorCode ierr;
598bf7f4e0aSPeter Brune 
599bf388a1fSBarry Smith   PetscFunctionBegin;
600f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
601bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
602bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
603*064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(Y,VEC_CLASSID,5);
604bf7f4e0aSPeter Brune 
605422a814eSBarry Smith   linesearch->result = SNES_LINESEARCH_SUCCEEDED;
606bf7f4e0aSPeter Brune 
607bf7f4e0aSPeter Brune   linesearch->vec_sol    = X;
608bf7f4e0aSPeter Brune   linesearch->vec_update = Y;
609bf7f4e0aSPeter Brune   linesearch->vec_func   = F;
610bf7f4e0aSPeter Brune 
611f1c6b773SPeter Brune   ierr = SNESLineSearchSetUp(linesearch);CHKERRQ(ierr);
612bf7f4e0aSPeter Brune 
613f5af7f23SKarl Rupp   if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */
614bf7f4e0aSPeter Brune 
615f5af7f23SKarl Rupp   if (fnorm) linesearch->fnorm = *fnorm;
616f5af7f23SKarl Rupp   else {
617bf7f4e0aSPeter Brune     ierr = VecNorm(F, NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
618bf7f4e0aSPeter Brune   }
619bf7f4e0aSPeter Brune 
62057a83faaSBarry Smith   ierr = PetscLogEventBegin(SNESLINESEARCH_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
621bf7f4e0aSPeter Brune 
622bf7f4e0aSPeter Brune   ierr = (*linesearch->ops->apply)(linesearch);CHKERRQ(ierr);
623bf7f4e0aSPeter Brune 
62457a83faaSBarry Smith   ierr = PetscLogEventEnd(SNESLINESEARCH_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
625bf7f4e0aSPeter Brune 
626f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
627bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
628bf7f4e0aSPeter Brune }
629bf7f4e0aSPeter Brune 
630f40b411bSPeter Brune /*@
631f1c6b773SPeter Brune    SNESLineSearchDestroy - Destroys the line search instance.
632f40b411bSPeter Brune 
633f1c6b773SPeter Brune    Collective on SNESLineSearch
634f40b411bSPeter Brune 
635f40b411bSPeter Brune    Input Parameters:
6368e557f58SPeter Brune .  linesearch - The linesearch context
637f40b411bSPeter Brune 
63884238204SBarry Smith    Level: developer
639f40b411bSPeter Brune 
64084238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchCreate(), SNESLineSearchReset(), SNESDestroy()
641f40b411bSPeter Brune @*/
642bf388a1fSBarry Smith PetscErrorCode SNESLineSearchDestroy(SNESLineSearch * linesearch)
643bf388a1fSBarry Smith {
644bf7f4e0aSPeter Brune   PetscErrorCode ierr;
645bf388a1fSBarry Smith 
646bf7f4e0aSPeter Brune   PetscFunctionBegin;
647bf7f4e0aSPeter Brune   if (!*linesearch) PetscFunctionReturn(0);
648f1c6b773SPeter Brune   PetscValidHeaderSpecific((*linesearch),SNESLINESEARCH_CLASSID,1);
6499e5d0892SLisandro Dalcin   if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = NULL; PetscFunctionReturn(0);}
650e04113cfSBarry Smith   ierr = PetscObjectSAWsViewOff((PetscObject)*linesearch);CHKERRQ(ierr);
65122d28d08SBarry Smith   ierr = SNESLineSearchReset(*linesearch);CHKERRQ(ierr);
652f5af7f23SKarl Rupp   if ((*linesearch)->ops->destroy) (*linesearch)->ops->destroy(*linesearch);
653bf7f4e0aSPeter Brune   ierr = PetscViewerDestroy(&(*linesearch)->monitor);CHKERRQ(ierr);
654dcf2fd19SBarry Smith   ierr = SNESLineSearchMonitorCancel((*linesearch));CHKERRQ(ierr);
655e7058c64SPeter Brune   ierr = PetscHeaderDestroy(linesearch);CHKERRQ(ierr);
656bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
657bf7f4e0aSPeter Brune }
658bf7f4e0aSPeter Brune 
659f40b411bSPeter Brune /*@
660dcf2fd19SBarry Smith    SNESLineSearchSetDefaultMonitor - Turns on/off printing useful information and debugging output about the line search.
661bf7f4e0aSPeter Brune 
662bf7f4e0aSPeter Brune    Input Parameters:
663dcf2fd19SBarry Smith +  linesearch - the linesearch object
664dcf2fd19SBarry Smith -  viewer - an ASCII PetscViewer or NULL to turn off monitor
665bf7f4e0aSPeter Brune 
666dcf2fd19SBarry Smith    Logically Collective on SNESLineSearch
667bf7f4e0aSPeter Brune 
668bf7f4e0aSPeter Brune    Options Database:
669dcf2fd19SBarry Smith .   -snes_linesearch_monitor [:filename] - enables the monitor
670bf7f4e0aSPeter Brune 
671bf7f4e0aSPeter Brune    Level: intermediate
672bf7f4e0aSPeter Brune 
673dcf2fd19SBarry Smith    Developer Note: This monitor is implemented differently than the other SNESLineSearchMonitors that are set with
674d12e167eSBarry Smith      SNESLineSearchMonitorSet() since it is called in many locations of the line search routines to display aspects of the
675d12e167eSBarry Smith      line search that are not visible to the other monitors.
676dcf2fd19SBarry Smith 
67784238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchGetDefaultMonitor(), PetscViewer, SNESLineSearchSetMonitor()
678bf7f4e0aSPeter Brune @*/
679dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchSetDefaultMonitor(SNESLineSearch linesearch, PetscViewer viewer)
680bf7f4e0aSPeter Brune {
681bf7f4e0aSPeter Brune   PetscErrorCode ierr;
682bf388a1fSBarry Smith 
683bf7f4e0aSPeter Brune   PetscFunctionBegin;
684dcf2fd19SBarry Smith   if (viewer) {ierr = PetscObjectReference((PetscObject)viewer);CHKERRQ(ierr);}
685bf7f4e0aSPeter Brune   ierr = PetscViewerDestroy(&linesearch->monitor);CHKERRQ(ierr);
686dcf2fd19SBarry Smith   linesearch->monitor = viewer;
687bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
688bf7f4e0aSPeter Brune }
689bf7f4e0aSPeter Brune 
690f40b411bSPeter Brune /*@
691dcf2fd19SBarry Smith    SNESLineSearchGetDefaultMonitor - Gets the PetscViewer instance for the line search monitor.
6926a388c36SPeter Brune 
693f190f2fcSBarry Smith    Input Parameter:
6948e557f58SPeter Brune .  linesearch - linesearch context
695f40b411bSPeter Brune 
696f190f2fcSBarry Smith    Output Parameter:
6978e557f58SPeter Brune .  monitor - monitor context
698f40b411bSPeter Brune 
699f40b411bSPeter Brune    Logically Collective on SNES
700f40b411bSPeter Brune 
701f40b411bSPeter Brune    Options Database Keys:
7028e557f58SPeter Brune .   -snes_linesearch_monitor - enables the monitor
703f40b411bSPeter Brune 
704f40b411bSPeter Brune    Level: intermediate
705f40b411bSPeter Brune 
70684238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchSetDefaultMonitor(), PetscViewer
707f40b411bSPeter Brune @*/
708dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchGetDefaultMonitor(SNESLineSearch linesearch, PetscViewer *monitor)
7096a388c36SPeter Brune {
7106a388c36SPeter Brune   PetscFunctionBegin;
711f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
7126a388c36SPeter Brune   *monitor = linesearch->monitor;
7136a388c36SPeter Brune   PetscFunctionReturn(0);
7146a388c36SPeter Brune }
7156a388c36SPeter Brune 
716dcf2fd19SBarry Smith /*@C
717dcf2fd19SBarry Smith    SNESLineSearchMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
718dcf2fd19SBarry Smith 
719dcf2fd19SBarry Smith    Collective on SNESLineSearch
720dcf2fd19SBarry Smith 
721dcf2fd19SBarry Smith    Input Parameters:
722dcf2fd19SBarry Smith +  ls - LineSearch object you wish to monitor
723dcf2fd19SBarry Smith .  name - the monitor type one is seeking
724dcf2fd19SBarry Smith .  help - message indicating what monitoring is done
725dcf2fd19SBarry Smith .  manual - manual page for the monitor
726dcf2fd19SBarry Smith .  monitor - the monitor function
727dcf2fd19SBarry 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
728dcf2fd19SBarry Smith 
729dcf2fd19SBarry Smith    Level: developer
730dcf2fd19SBarry Smith 
731dcf2fd19SBarry Smith .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
732dcf2fd19SBarry Smith           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
733dcf2fd19SBarry Smith           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
734dcf2fd19SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
735dcf2fd19SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
736dcf2fd19SBarry Smith           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
737dcf2fd19SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
738dcf2fd19SBarry Smith @*/
739d12e167eSBarry Smith PetscErrorCode  SNESLineSearchMonitorSetFromOptions(SNESLineSearch ls,const char name[],const char help[], const char manual[],PetscErrorCode (*monitor)(SNESLineSearch,PetscViewerAndFormat*),PetscErrorCode (*monitorsetup)(SNESLineSearch,PetscViewerAndFormat*))
740dcf2fd19SBarry Smith {
741dcf2fd19SBarry Smith   PetscErrorCode    ierr;
742dcf2fd19SBarry Smith   PetscViewer       viewer;
743dcf2fd19SBarry Smith   PetscViewerFormat format;
744dcf2fd19SBarry Smith   PetscBool         flg;
745dcf2fd19SBarry Smith 
746dcf2fd19SBarry Smith   PetscFunctionBegin;
74716413a6aSBarry Smith   ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)ls),((PetscObject) ls)->options,((PetscObject)ls)->prefix,name,&viewer,&format,&flg);CHKERRQ(ierr);
748dcf2fd19SBarry Smith   if (flg) {
749d12e167eSBarry Smith     PetscViewerAndFormat *vf;
750d12e167eSBarry Smith     ierr = PetscViewerAndFormatCreate(viewer,format,&vf);CHKERRQ(ierr);
751d12e167eSBarry Smith     ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr);
752dcf2fd19SBarry Smith     if (monitorsetup) {
753d12e167eSBarry Smith       ierr = (*monitorsetup)(ls,vf);CHKERRQ(ierr);
754dcf2fd19SBarry Smith     }
755d12e167eSBarry Smith     ierr = SNESLineSearchMonitorSet(ls,(PetscErrorCode (*)(SNESLineSearch,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);CHKERRQ(ierr);
756dcf2fd19SBarry Smith   }
757dcf2fd19SBarry Smith   PetscFunctionReturn(0);
758dcf2fd19SBarry Smith }
759dcf2fd19SBarry Smith 
760f40b411bSPeter Brune /*@
761f1c6b773SPeter Brune    SNESLineSearchSetFromOptions - Sets options for the line search
762f40b411bSPeter Brune 
763f40b411bSPeter Brune    Input Parameters:
7648e557f58SPeter Brune .  linesearch - linesearch context
765f40b411bSPeter Brune 
766f40b411bSPeter Brune    Options Database Keys:
767d4c6564cSPatrick Farrell + -snes_linesearch_type <type> - basic, bt, l2, cp, nleqerr, shell
7683c7d6663SPeter Brune . -snes_linesearch_order <order> - 1, 2, 3.  Most types only support certain orders (bt supports 2 or 3)
7691fe24845SBarry Smith . -snes_linesearch_norms   - Turn on/off the linesearch norms for the basic linesearch typem (SNESLineSearchSetComputeNorms())
77071eef1aeSPeter Brune . -snes_linesearch_minlambda - The minimum step length
7711a9b3a06SPeter Brune . -snes_linesearch_maxstep - The maximum step size
7721a9b3a06SPeter Brune . -snes_linesearch_rtol - Relative tolerance for iterative line searches
7731a9b3a06SPeter Brune . -snes_linesearch_atol - Absolute tolerance for iterative line searches
7741a9b3a06SPeter Brune . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches
7751a9b3a06SPeter Brune . -snes_linesearch_max_it - The number of iterations for iterative line searches
776dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches
777dcf2fd19SBarry Smith . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine
7788e557f58SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter
779cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess.
7801a9b3a06SPeter Brune . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method
781d8d34be6SBarry Smith - -snes_linesearch_precheck_picard_angle - Angle used in Picard precheck method
782f40b411bSPeter Brune 
783f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
784f40b411bSPeter Brune 
785f40b411bSPeter Brune    Level: intermediate
786f40b411bSPeter Brune 
7871fe24845SBarry Smith .seealso: SNESLineSearchCreate(), SNESLineSearchSetOrder(), SNESLineSearchSetType(), SNESLineSearchSetTolerances(), SNESLineSearchSetDamping(), SNESLineSearchPreCheckPicard(),
7881fe24845SBarry Smith           SNESLineSearchType, SNESLineSearchSetComputeNorms()
789f40b411bSPeter Brune @*/
790bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch)
791bf388a1fSBarry Smith {
792bf7f4e0aSPeter Brune   PetscErrorCode    ierr;
7931a4f838cSPeter Brune   const char        *deft = SNESLINESEARCHBASIC;
794bf7f4e0aSPeter Brune   char              type[256];
795bf7f4e0aSPeter Brune   PetscBool         flg, set;
796dcf2fd19SBarry Smith   PetscViewer       viewer;
797bf388a1fSBarry Smith 
798bf7f4e0aSPeter Brune   PetscFunctionBegin;
7990f51fdf8SToby Isaac   ierr = SNESLineSearchRegisterAll();CHKERRQ(ierr);
800bf7f4e0aSPeter Brune 
801bf7f4e0aSPeter Brune   ierr = PetscObjectOptionsBegin((PetscObject)linesearch);CHKERRQ(ierr);
802f5af7f23SKarl Rupp   if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name;
803a264d7a6SBarry Smith   ierr = PetscOptionsFList("-snes_linesearch_type","Linesearch type","SNESLineSearchSetType",SNESLineSearchList,deft,type,256,&flg);CHKERRQ(ierr);
804bf7f4e0aSPeter Brune   if (flg) {
805f1c6b773SPeter Brune     ierr = SNESLineSearchSetType(linesearch,type);CHKERRQ(ierr);
806bf7f4e0aSPeter Brune   } else if (!((PetscObject)linesearch)->type_name) {
807f1c6b773SPeter Brune     ierr = SNESLineSearchSetType(linesearch,deft);CHKERRQ(ierr);
808bf7f4e0aSPeter Brune   }
809bf7f4e0aSPeter Brune 
81016413a6aSBarry Smith   ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)linesearch),((PetscObject) linesearch)->options,((PetscObject)linesearch)->prefix,"-snes_linesearch_monitor",&viewer,NULL,&set);CHKERRQ(ierr);
811dcf2fd19SBarry Smith   if (set) {
812dcf2fd19SBarry Smith     ierr = SNESLineSearchSetDefaultMonitor(linesearch,viewer);CHKERRQ(ierr);
813dcf2fd19SBarry Smith     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
814dcf2fd19SBarry Smith   }
815dcf2fd19SBarry Smith   ierr = SNESLineSearchMonitorSetFromOptions(linesearch,"-snes_linesearch_monitor_solution_update","View correction at each iteration","SNESLineSearchMonitorSolutionUpdate",SNESLineSearchMonitorSolutionUpdate,NULL);CHKERRQ(ierr);
816bf7f4e0aSPeter Brune 
8171a9b3a06SPeter Brune   /* tolerances */
81894ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_minlambda","Minimum step length","SNESLineSearchSetTolerances",linesearch->steptol,&linesearch->steptol,NULL);CHKERRQ(ierr);
81994ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_maxstep","Maximum step size","SNESLineSearchSetTolerances",linesearch->maxstep,&linesearch->maxstep,NULL);CHKERRQ(ierr);
82094ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_rtol","Relative tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->rtol,&linesearch->rtol,NULL);CHKERRQ(ierr);
82194ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_atol","Absolute tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->atol,&linesearch->atol,NULL);CHKERRQ(ierr);
82294ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_ltol","Change in lambda tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->ltol,&linesearch->ltol,NULL);CHKERRQ(ierr);
82394ae4db5SBarry Smith   ierr = PetscOptionsInt("-snes_linesearch_max_it","Maximum iterations for iterative line searches","SNESLineSearchSetTolerances",linesearch->max_its,&linesearch->max_its,NULL);CHKERRQ(ierr);
8247a35526eSPeter Brune 
8251a9b3a06SPeter Brune   /* damping parameters */
82694ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_damping","Line search damping and initial step guess","SNESLineSearchSetDamping",linesearch->damping,&linesearch->damping,NULL);CHKERRQ(ierr);
8271a9b3a06SPeter Brune 
82894ae4db5SBarry Smith   ierr = PetscOptionsBool("-snes_linesearch_keeplambda","Use previous lambda as damping","SNESLineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,NULL);CHKERRQ(ierr);
8291a9b3a06SPeter Brune 
8301a9b3a06SPeter Brune   /* precheck */
8317a35526eSPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_precheck_picard","Use a correction that sometimes improves convergence of Picard iteration","SNESLineSearchPreCheckPicard",flg,&flg,&set);CHKERRQ(ierr);
8327a35526eSPeter Brune   if (set) {
8337a35526eSPeter Brune     if (flg) {
8347a35526eSPeter Brune       linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */
835f5af7f23SKarl Rupp 
8367a35526eSPeter Brune       ierr = PetscOptionsReal("-snes_linesearch_precheck_picard_angle","Maximum angle at which to activate the correction",
8370298fd71SBarry Smith                               "none",linesearch->precheck_picard_angle,&linesearch->precheck_picard_angle,NULL);CHKERRQ(ierr);
8387a35526eSPeter Brune       ierr = SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&linesearch->precheck_picard_angle);CHKERRQ(ierr);
8397a35526eSPeter Brune     } else {
8400298fd71SBarry Smith       ierr = SNESLineSearchSetPreCheck(linesearch,NULL,NULL);CHKERRQ(ierr);
8417a35526eSPeter Brune     }
8427a35526eSPeter Brune   }
84394ae4db5SBarry Smith   ierr = PetscOptionsInt("-snes_linesearch_order","Order of approximation used in the line search","SNESLineSearchSetOrder",linesearch->order,&linesearch->order,NULL);CHKERRQ(ierr);
84494ae4db5SBarry Smith   ierr = PetscOptionsBool("-snes_linesearch_norms","Compute final norms in line search","SNESLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,NULL);CHKERRQ(ierr);
8457a35526eSPeter Brune 
8465a9b6599SPeter Brune   if (linesearch->ops->setfromoptions) {
847e55864a3SBarry Smith     (*linesearch->ops->setfromoptions)(PetscOptionsObject,linesearch);CHKERRQ(ierr);
8485a9b6599SPeter Brune   }
8495a9b6599SPeter Brune 
8500633abcbSJed Brown   ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)linesearch);CHKERRQ(ierr);
851bf7f4e0aSPeter Brune   ierr = PetscOptionsEnd();CHKERRQ(ierr);
852bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
853bf7f4e0aSPeter Brune }
854bf7f4e0aSPeter Brune 
855f40b411bSPeter Brune /*@
856f190f2fcSBarry Smith    SNESLineSearchView - Prints useful information about the line search
857f40b411bSPeter Brune 
858f40b411bSPeter Brune    Input Parameters:
8598e557f58SPeter Brune .  linesearch - linesearch context
860f40b411bSPeter Brune 
861f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
862f40b411bSPeter Brune 
863f40b411bSPeter Brune    Level: intermediate
864f40b411bSPeter Brune 
865dcf2fd19SBarry Smith .seealso: SNESLineSearchCreate()
866f40b411bSPeter Brune @*/
867bf388a1fSBarry Smith PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer)
868bf388a1fSBarry Smith {
8697f1410a3SPeter Brune   PetscErrorCode ierr;
8707f1410a3SPeter Brune   PetscBool      iascii;
871bf388a1fSBarry Smith 
872bf7f4e0aSPeter Brune   PetscFunctionBegin;
8737f1410a3SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
8747f1410a3SPeter Brune   if (!viewer) {
875ce94432eSBarry Smith     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch),&viewer);CHKERRQ(ierr);
8767f1410a3SPeter Brune   }
8777f1410a3SPeter Brune   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
8787f1410a3SPeter Brune   PetscCheckSameComm(linesearch,1,viewer,2);
879f40b411bSPeter Brune 
880251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
8817f1410a3SPeter Brune   if (iascii) {
882dae58748SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)linesearch,viewer);CHKERRQ(ierr);
8837f1410a3SPeter Brune     if (linesearch->ops->view) {
8847f1410a3SPeter Brune       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
8857f1410a3SPeter Brune       ierr = (*linesearch->ops->view)(linesearch,viewer);CHKERRQ(ierr);
8867f1410a3SPeter Brune       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
8877f1410a3SPeter Brune     }
888c69d1a72SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep,(double)linesearch->steptol);CHKERRQ(ierr);
889c69d1a72SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol,(double)linesearch->atol,(double)linesearch->ltol);CHKERRQ(ierr);
8907f1410a3SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D\n", linesearch->max_its);CHKERRQ(ierr);
8916b2b7091SBarry Smith     if (linesearch->ops->precheck) {
8926b2b7091SBarry Smith       if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) {
8937f1410a3SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  using precheck step to speed up Picard convergence\n", linesearch->max_its);CHKERRQ(ierr);
8947f1410a3SPeter Brune       } else {
8957f1410a3SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  using user-defined precheck step\n", linesearch->max_its);CHKERRQ(ierr);
8967f1410a3SPeter Brune       }
8977f1410a3SPeter Brune     }
8986b2b7091SBarry Smith     if (linesearch->ops->postcheck) {
8997f1410a3SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer,"  using user-defined postcheck step\n", linesearch->max_its);CHKERRQ(ierr);
9007f1410a3SPeter Brune     }
9017f1410a3SPeter Brune   }
902bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
903bf7f4e0aSPeter Brune }
904bf7f4e0aSPeter Brune 
905ea5d4fccSPeter Brune /*@C
906a80ff896SJed Brown    SNESLineSearchGetType - Gets the linesearch type
907a80ff896SJed Brown 
908a80ff896SJed Brown    Logically Collective on SNESLineSearch
909a80ff896SJed Brown 
910a80ff896SJed Brown    Input Parameters:
911a80ff896SJed Brown .  linesearch - linesearch context
912a80ff896SJed Brown 
913a80ff896SJed Brown    Output Parameters:
914a80ff896SJed Brown -  type - The type of line search, or NULL if not set
915a80ff896SJed Brown 
916a80ff896SJed Brown    Level: intermediate
917a80ff896SJed Brown 
918a80ff896SJed Brown .seealso: SNESLineSearchCreate(), SNESLineSearchType, SNESLineSearchSetFromOptions(), SNESLineSearchSetType()
919a80ff896SJed Brown @*/
920a80ff896SJed Brown PetscErrorCode SNESLineSearchGetType(SNESLineSearch linesearch, SNESLineSearchType *type)
921a80ff896SJed Brown {
922a80ff896SJed Brown   PetscFunctionBegin;
923a80ff896SJed Brown   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
924a80ff896SJed Brown   PetscValidCharPointer(type,2);
925a80ff896SJed Brown   *type = ((PetscObject)linesearch)->type_name;
926a80ff896SJed Brown   PetscFunctionReturn(0);
927a80ff896SJed Brown }
928a80ff896SJed Brown 
929a80ff896SJed Brown /*@C
930f1c6b773SPeter Brune    SNESLineSearchSetType - Sets the linesearch type
931f40b411bSPeter Brune 
932f190f2fcSBarry Smith    Logically Collective on SNESLineSearch
933f190f2fcSBarry Smith 
934f40b411bSPeter Brune    Input Parameters:
9358e557f58SPeter Brune +  linesearch - linesearch context
936f40b411bSPeter Brune -  type - The type of line search to be used
937f40b411bSPeter Brune 
938cd7522eaSPeter Brune    Available Types:
9391fe24845SBarry Smith +  SNESLINESEARCHBASIC - Simple damping line search, defaults to using the full Newton step
9401fe24845SBarry Smith .  SNESLINESEARCHBT - Backtracking line search over the L2 norm of the function
9411fe24845SBarry Smith .  SNESLINESEARCHL2 - Secant line search over the L2 norm of the function
9421fe24845SBarry Smith .  SNESLINESEARCHCP - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x)
9431fe24845SBarry Smith .  SNESLINESEARCHNLEQERR - Affine-covariant error-oriented linesearch
9441fe24845SBarry Smith -  SNESLINESEARCHSHELL - User provided SNESLineSearch implementation
9451fe24845SBarry Smith 
9461fe24845SBarry Smith    Options Database:
9471fe24845SBarry Smith .  -snes_linesearch_type <type> - basic, bt, l2, cp, nleqerr, shell
948cd7522eaSPeter Brune 
949f40b411bSPeter Brune    Level: intermediate
950f40b411bSPeter Brune 
951a80ff896SJed Brown .seealso: SNESLineSearchCreate(), SNESLineSearchType, SNESLineSearchSetFromOptions(), SNESLineSearchGetType()
952f40b411bSPeter Brune @*/
95319fd82e9SBarry Smith PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type)
954bf7f4e0aSPeter Brune {
955f1c6b773SPeter Brune   PetscErrorCode ierr,(*r)(SNESLineSearch);
956bf7f4e0aSPeter Brune   PetscBool      match;
957bf7f4e0aSPeter Brune 
958bf7f4e0aSPeter Brune   PetscFunctionBegin;
959f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
960bf7f4e0aSPeter Brune   PetscValidCharPointer(type,2);
961bf7f4e0aSPeter Brune 
962251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)linesearch,type,&match);CHKERRQ(ierr);
963bf7f4e0aSPeter Brune   if (match) PetscFunctionReturn(0);
964bf7f4e0aSPeter Brune 
9651c9cd337SJed Brown   ierr = PetscFunctionListFind(SNESLineSearchList,type,&r);CHKERRQ(ierr);
966bf7f4e0aSPeter Brune   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type);
967bf7f4e0aSPeter Brune   /* Destroy the previous private linesearch context */
968bf7f4e0aSPeter Brune   if (linesearch->ops->destroy) {
969bf7f4e0aSPeter Brune     ierr = (*(linesearch)->ops->destroy)(linesearch);CHKERRQ(ierr);
9700298fd71SBarry Smith     linesearch->ops->destroy = NULL;
971bf7f4e0aSPeter Brune   }
972f1c6b773SPeter Brune   /* Reinitialize function pointers in SNESLineSearchOps structure */
9739e5d0892SLisandro Dalcin   linesearch->ops->apply          = NULL;
9749e5d0892SLisandro Dalcin   linesearch->ops->view           = NULL;
9759e5d0892SLisandro Dalcin   linesearch->ops->setfromoptions = NULL;
9769e5d0892SLisandro Dalcin   linesearch->ops->destroy        = NULL;
977bf7f4e0aSPeter Brune 
978bf7f4e0aSPeter Brune   ierr = PetscObjectChangeTypeName((PetscObject)linesearch,type);CHKERRQ(ierr);
979bf7f4e0aSPeter Brune   ierr = (*r)(linesearch);CHKERRQ(ierr);
980bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
981bf7f4e0aSPeter Brune }
982bf7f4e0aSPeter Brune 
983f40b411bSPeter Brune /*@
98478bcb3b5SPeter Brune    SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation.
985f40b411bSPeter Brune 
986f40b411bSPeter Brune    Input Parameters:
9878e557f58SPeter Brune +  linesearch - linesearch context
988f40b411bSPeter Brune -  snes - The snes instance
989f40b411bSPeter Brune 
99078bcb3b5SPeter Brune    Level: developer
99178bcb3b5SPeter Brune 
99278bcb3b5SPeter Brune    Notes:
993f190f2fcSBarry Smith    This happens automatically when the line search is obtained/created with
9947601faf0SJed Brown    SNESGetLineSearch().  This routine is therefore mainly called within SNES
99578bcb3b5SPeter Brune    implementations.
996f40b411bSPeter Brune 
9978141a3b9SPeter Brune    Level: developer
998f40b411bSPeter Brune 
999cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
1000f40b411bSPeter Brune @*/
1001bf388a1fSBarry Smith PetscErrorCode  SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes)
1002bf388a1fSBarry Smith {
1003bf7f4e0aSPeter Brune   PetscFunctionBegin;
1004f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1005bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,2);
1006bf7f4e0aSPeter Brune   linesearch->snes = snes;
1007bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1008bf7f4e0aSPeter Brune }
1009bf7f4e0aSPeter Brune 
1010f40b411bSPeter Brune /*@
10118141a3b9SPeter Brune    SNESLineSearchGetSNES - Gets the SNES instance associated with the line search.
10128141a3b9SPeter Brune    Having an associated SNES is necessary because most line search implementations must be able to
10138141a3b9SPeter Brune    evaluate the function using SNESComputeFunction() for the associated SNES.  This routine
10148141a3b9SPeter Brune    is used in the line search implementations when one must get this associated SNES instance.
1015f40b411bSPeter Brune 
1016f40b411bSPeter Brune    Input Parameters:
10178e557f58SPeter Brune .  linesearch - linesearch context
1018f40b411bSPeter Brune 
1019f40b411bSPeter Brune    Output Parameters:
1020f40b411bSPeter Brune .  snes - The snes instance
1021f40b411bSPeter Brune 
10228141a3b9SPeter Brune    Level: developer
1023f40b411bSPeter Brune 
1024cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
1025f40b411bSPeter Brune @*/
1026bf388a1fSBarry Smith PetscErrorCode  SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes)
1027bf388a1fSBarry Smith {
1028bf7f4e0aSPeter Brune   PetscFunctionBegin;
1029f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
10306a388c36SPeter Brune   PetscValidPointer(snes, 2);
1031bf7f4e0aSPeter Brune   *snes = linesearch->snes;
1032bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1033bf7f4e0aSPeter Brune }
1034bf7f4e0aSPeter Brune 
1035f40b411bSPeter Brune /*@
1036f1c6b773SPeter Brune    SNESLineSearchGetLambda - Gets the last linesearch steplength discovered.
1037f40b411bSPeter Brune 
1038f40b411bSPeter Brune    Input Parameters:
10398e557f58SPeter Brune .  linesearch - linesearch context
1040f40b411bSPeter Brune 
1041f40b411bSPeter Brune    Output Parameters:
1042cd7522eaSPeter Brune .  lambda - The last steplength computed during SNESLineSearchApply()
1043f40b411bSPeter Brune 
104478bcb3b5SPeter Brune    Level: advanced
104578bcb3b5SPeter Brune 
10468e557f58SPeter Brune    Notes:
10478e557f58SPeter Brune    This is useful in methods where the solver is ill-scaled and
104878bcb3b5SPeter Brune    requires some adaptive notion of the difference in scale between the
104978bcb3b5SPeter Brune    solution and the function.  For instance, SNESQN may be scaled by the
105078bcb3b5SPeter Brune    line search lambda using the argument -snes_qn_scaling ls.
105178bcb3b5SPeter Brune 
1052cd7522eaSPeter Brune .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetDamping(), SNESLineSearchApply()
1053f40b411bSPeter Brune @*/
1054f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda)
10556a388c36SPeter Brune {
10566a388c36SPeter Brune   PetscFunctionBegin;
1057f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1058534a8f05SLisandro Dalcin   PetscValidRealPointer(lambda, 2);
10596a388c36SPeter Brune   *lambda = linesearch->lambda;
10606a388c36SPeter Brune   PetscFunctionReturn(0);
10616a388c36SPeter Brune }
10626a388c36SPeter Brune 
1063f40b411bSPeter Brune /*@
1064f1c6b773SPeter Brune    SNESLineSearchSetLambda - Sets the linesearch steplength.
1065f40b411bSPeter Brune 
1066f40b411bSPeter Brune    Input Parameters:
10678e557f58SPeter Brune +  linesearch - linesearch context
1068f40b411bSPeter Brune -  lambda - The last steplength.
1069f40b411bSPeter Brune 
1070cd7522eaSPeter Brune    Notes:
1071f190f2fcSBarry Smith    This routine is typically used within implementations of SNESLineSearchApply()
1072cd7522eaSPeter Brune    to set the final steplength.  This routine (and SNESLineSearchGetLambda()) were
1073cd7522eaSPeter Brune    added in order to facilitate Quasi-Newton methods that use the previous steplength
1074cd7522eaSPeter Brune    as an inner scaling parameter.
1075cd7522eaSPeter Brune 
107678bcb3b5SPeter Brune    Level: advanced
1077f40b411bSPeter Brune 
1078f1c6b773SPeter Brune .seealso: SNESLineSearchGetLambda()
1079f40b411bSPeter Brune @*/
1080f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
10816a388c36SPeter Brune {
10826a388c36SPeter Brune   PetscFunctionBegin;
1083f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
10846a388c36SPeter Brune   linesearch->lambda = lambda;
10856a388c36SPeter Brune   PetscFunctionReturn(0);
10866a388c36SPeter Brune }
10876a388c36SPeter Brune 
1088f40b411bSPeter Brune /*@
10893c7d6663SPeter Brune    SNESLineSearchGetTolerances - Gets the tolerances for the linesearch.  These include
109078bcb3b5SPeter Brune    tolerances for the relative and absolute change in the function norm, the change
109178bcb3b5SPeter Brune    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
109278bcb3b5SPeter Brune    and the maximum number of iterations the line search procedure may take.
1093f40b411bSPeter Brune 
1094f40b411bSPeter Brune    Input Parameters:
10958e557f58SPeter Brune .  linesearch - linesearch context
1096f40b411bSPeter Brune 
1097f40b411bSPeter Brune    Output Parameters:
1098516fe3c3SPeter Brune +  steptol - The minimum steplength
10996cc8e53bSPeter Brune .  maxstep - The maximum steplength
1100516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
1101516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
1102516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
1103516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
1104f40b411bSPeter Brune 
110578bcb3b5SPeter Brune    Level: intermediate
110678bcb3b5SPeter Brune 
110778bcb3b5SPeter Brune    Notes:
110878bcb3b5SPeter Brune    Different line searches may implement these parameters slightly differently as
11093c7d6663SPeter Brune    the type requires.
1110516fe3c3SPeter Brune 
1111f1c6b773SPeter Brune .seealso: SNESLineSearchSetTolerances()
1112f40b411bSPeter Brune @*/
1113f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
11146a388c36SPeter Brune {
11156a388c36SPeter Brune   PetscFunctionBegin;
1116f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1117516fe3c3SPeter Brune   if (steptol) {
1118534a8f05SLisandro Dalcin     PetscValidRealPointer(steptol, 2);
11196a388c36SPeter Brune     *steptol = linesearch->steptol;
1120516fe3c3SPeter Brune   }
1121516fe3c3SPeter Brune   if (maxstep) {
1122534a8f05SLisandro Dalcin     PetscValidRealPointer(maxstep, 3);
1123516fe3c3SPeter Brune     *maxstep = linesearch->maxstep;
1124516fe3c3SPeter Brune   }
1125516fe3c3SPeter Brune   if (rtol) {
1126534a8f05SLisandro Dalcin     PetscValidRealPointer(rtol, 4);
1127516fe3c3SPeter Brune     *rtol = linesearch->rtol;
1128516fe3c3SPeter Brune   }
1129516fe3c3SPeter Brune   if (atol) {
1130534a8f05SLisandro Dalcin     PetscValidRealPointer(atol, 5);
1131516fe3c3SPeter Brune     *atol = linesearch->atol;
1132516fe3c3SPeter Brune   }
1133516fe3c3SPeter Brune   if (ltol) {
1134534a8f05SLisandro Dalcin     PetscValidRealPointer(ltol, 6);
1135516fe3c3SPeter Brune     *ltol = linesearch->ltol;
1136516fe3c3SPeter Brune   }
1137516fe3c3SPeter Brune   if (max_its) {
1138534a8f05SLisandro Dalcin     PetscValidIntPointer(max_its, 7);
1139516fe3c3SPeter Brune     *max_its = linesearch->max_its;
1140516fe3c3SPeter Brune   }
11416a388c36SPeter Brune   PetscFunctionReturn(0);
11426a388c36SPeter Brune }
11436a388c36SPeter Brune 
1144f40b411bSPeter Brune /*@
11453c7d6663SPeter Brune    SNESLineSearchSetTolerances -  Gets the tolerances for the linesearch.  These include
114678bcb3b5SPeter Brune    tolerances for the relative and absolute change in the function norm, the change
114778bcb3b5SPeter Brune    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
114878bcb3b5SPeter Brune    and the maximum number of iterations the line search procedure may take.
1149f40b411bSPeter Brune 
1150f40b411bSPeter Brune    Input Parameters:
11518e557f58SPeter Brune +  linesearch - linesearch context
1152516fe3c3SPeter Brune .  steptol - The minimum steplength
11536cc8e53bSPeter Brune .  maxstep - The maximum steplength
1154516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
1155516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
1156516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
1157516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
1158f40b411bSPeter Brune 
115978bcb3b5SPeter Brune    Notes:
11603c7d6663SPeter Brune    The user may choose to not set any of the tolerances using PETSC_DEFAULT in
116178bcb3b5SPeter Brune    place of an argument.
1162f40b411bSPeter Brune 
116378bcb3b5SPeter Brune    Level: intermediate
1164516fe3c3SPeter Brune 
1165f1c6b773SPeter Brune .seealso: SNESLineSearchGetTolerances()
1166f40b411bSPeter Brune @*/
1167f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its)
11686a388c36SPeter Brune {
11696a388c36SPeter Brune   PetscFunctionBegin;
1170f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1171d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,steptol,2);
1172d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,maxstep,3);
1173d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,rtol,4);
1174d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,atol,5);
1175d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,ltol,6);
1176d3952184SSatish Balay   PetscValidLogicalCollectiveInt(linesearch,max_its,7);
1177d3952184SSatish Balay 
1178d3952184SSatish Balay   if (steptol!= PETSC_DEFAULT) {
1179ce94432eSBarry Smith     if (steptol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Minimum step length %14.12e must be non-negative",(double)steptol);
11806a388c36SPeter Brune     linesearch->steptol = steptol;
1181d3952184SSatish Balay   }
1182d3952184SSatish Balay 
1183d3952184SSatish Balay   if (maxstep!= PETSC_DEFAULT) {
1184ce94432eSBarry Smith     if (maxstep < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum step length %14.12e must be non-negative",(double)maxstep);
1185516fe3c3SPeter Brune     linesearch->maxstep = maxstep;
1186d3952184SSatish Balay   }
1187d3952184SSatish Balay 
1188d3952184SSatish Balay   if (rtol != PETSC_DEFAULT) {
1189ce94432eSBarry 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);
1190516fe3c3SPeter Brune     linesearch->rtol = rtol;
1191d3952184SSatish Balay   }
1192d3952184SSatish Balay 
1193d3952184SSatish Balay   if (atol != PETSC_DEFAULT) {
1194ce94432eSBarry Smith     if (atol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %14.12e must be non-negative",(double)atol);
1195516fe3c3SPeter Brune     linesearch->atol = atol;
1196d3952184SSatish Balay   }
1197d3952184SSatish Balay 
1198d3952184SSatish Balay   if (ltol != PETSC_DEFAULT) {
11992479783cSJose E. Roman     if (ltol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Lambda tolerance %14.12e must be non-negative",(double)ltol);
1200516fe3c3SPeter Brune     linesearch->ltol = ltol;
1201d3952184SSatish Balay   }
1202d3952184SSatish Balay 
1203d3952184SSatish Balay   if (max_its != PETSC_DEFAULT) {
1204ce94432eSBarry Smith     if (max_its < 0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",max_its);
1205516fe3c3SPeter Brune     linesearch->max_its = max_its;
1206d3952184SSatish Balay   }
12076a388c36SPeter Brune   PetscFunctionReturn(0);
12086a388c36SPeter Brune }
12096a388c36SPeter Brune 
1210f40b411bSPeter Brune /*@
1211f1c6b773SPeter Brune    SNESLineSearchGetDamping - Gets the line search damping parameter.
1212f40b411bSPeter Brune 
1213f40b411bSPeter Brune    Input Parameters:
12148e557f58SPeter Brune .  linesearch - linesearch context
1215f40b411bSPeter Brune 
1216f40b411bSPeter Brune    Output Parameters:
12178e557f58SPeter Brune .  damping - The damping parameter
1218f40b411bSPeter Brune 
121978bcb3b5SPeter Brune    Level: advanced
1220f40b411bSPeter Brune 
122178bcb3b5SPeter Brune .seealso: SNESLineSearchGetStepTolerance(), SNESQN
1222f40b411bSPeter Brune @*/
1223f40b411bSPeter Brune 
1224f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping)
12256a388c36SPeter Brune {
12266a388c36SPeter Brune   PetscFunctionBegin;
1227f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1228534a8f05SLisandro Dalcin   PetscValidRealPointer(damping, 2);
12296a388c36SPeter Brune   *damping = linesearch->damping;
12306a388c36SPeter Brune   PetscFunctionReturn(0);
12316a388c36SPeter Brune }
12326a388c36SPeter Brune 
1233f40b411bSPeter Brune /*@
1234fd292e60Sprj-    SNESLineSearchSetDamping - Sets the line search damping parameter.
1235f40b411bSPeter Brune 
1236f40b411bSPeter Brune    Input Parameters:
123703fd4120SBarry Smith +  linesearch - linesearch context
123803fd4120SBarry Smith -  damping - The damping parameter
1239f40b411bSPeter Brune 
124003fd4120SBarry Smith    Options Database:
124103fd4120SBarry Smith .   -snes_linesearch_damping
1242f40b411bSPeter Brune    Level: intermediate
1243f40b411bSPeter Brune 
1244cd7522eaSPeter Brune    Notes:
1245cd7522eaSPeter Brune    The basic line search merely takes the update step scaled by the damping parameter.
1246cd7522eaSPeter Brune    The use of the damping parameter in the l2 and cp line searches is much more subtle;
124778bcb3b5SPeter Brune    it is used as a starting point in calculating the secant step. However, the eventual
1248cd7522eaSPeter Brune    step may be of greater length than the damping parameter.  In the bt line search it is
1249cd7522eaSPeter Brune    used as the maximum possible step length, as the bt line search only backtracks.
1250cd7522eaSPeter Brune 
1251f1c6b773SPeter Brune .seealso: SNESLineSearchGetDamping()
1252f40b411bSPeter Brune @*/
1253f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping)
12546a388c36SPeter Brune {
12556a388c36SPeter Brune   PetscFunctionBegin;
1256f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
12576a388c36SPeter Brune   linesearch->damping = damping;
12586a388c36SPeter Brune   PetscFunctionReturn(0);
12596a388c36SPeter Brune }
12606a388c36SPeter Brune 
126159405d5eSPeter Brune /*@
126259405d5eSPeter Brune    SNESLineSearchGetOrder - Gets the line search approximation order.
126359405d5eSPeter Brune 
126459405d5eSPeter Brune    Input Parameters:
126578bcb3b5SPeter Brune .  linesearch - linesearch context
126659405d5eSPeter Brune 
126759405d5eSPeter Brune    Output Parameters:
12688e557f58SPeter Brune .  order - The order
126959405d5eSPeter Brune 
127078bcb3b5SPeter Brune    Possible Values for order:
12713c7d6663SPeter Brune +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
12723c7d6663SPeter Brune .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
12733c7d6663SPeter Brune -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
127478bcb3b5SPeter Brune 
127559405d5eSPeter Brune    Level: intermediate
127659405d5eSPeter Brune 
127759405d5eSPeter Brune .seealso: SNESLineSearchSetOrder()
127859405d5eSPeter Brune @*/
127959405d5eSPeter Brune 
1280b000cd8dSPeter Brune PetscErrorCode  SNESLineSearchGetOrder(SNESLineSearch linesearch,PetscInt *order)
128159405d5eSPeter Brune {
128259405d5eSPeter Brune   PetscFunctionBegin;
128359405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1284534a8f05SLisandro Dalcin   PetscValidIntPointer(order, 2);
128559405d5eSPeter Brune   *order = linesearch->order;
128659405d5eSPeter Brune   PetscFunctionReturn(0);
128759405d5eSPeter Brune }
128859405d5eSPeter Brune 
128959405d5eSPeter Brune /*@
12901f8196a2SJed Brown    SNESLineSearchSetOrder - Sets the maximum order of the polynomial fit used in the line search
129159405d5eSPeter Brune 
129259405d5eSPeter Brune    Input Parameters:
1293a2b725a8SWilliam Gropp +  linesearch - linesearch context
1294a2b725a8SWilliam Gropp -  order - The damping parameter
129559405d5eSPeter Brune 
129659405d5eSPeter Brune    Level: intermediate
129759405d5eSPeter Brune 
129878bcb3b5SPeter Brune    Possible Values for order:
12993c7d6663SPeter Brune +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
13003c7d6663SPeter Brune .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
13013c7d6663SPeter Brune -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
130278bcb3b5SPeter Brune 
130359405d5eSPeter Brune    Notes:
130459405d5eSPeter Brune    Variable orders are supported by the following line searches:
130578bcb3b5SPeter Brune +  bt - cubic and quadratic
130678bcb3b5SPeter Brune -  cp - linear and quadratic
130759405d5eSPeter Brune 
13081f8196a2SJed Brown .seealso: SNESLineSearchGetOrder(), SNESLineSearchSetDamping()
130959405d5eSPeter Brune @*/
1310b000cd8dSPeter Brune PetscErrorCode  SNESLineSearchSetOrder(SNESLineSearch linesearch,PetscInt order)
131159405d5eSPeter Brune {
131259405d5eSPeter Brune   PetscFunctionBegin;
131359405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
131459405d5eSPeter Brune   linesearch->order = order;
131559405d5eSPeter Brune   PetscFunctionReturn(0);
131659405d5eSPeter Brune }
131759405d5eSPeter Brune 
1318f40b411bSPeter Brune /*@
1319f1c6b773SPeter Brune    SNESLineSearchGetNorms - Gets the norms for for X, Y, and F.
1320f40b411bSPeter Brune 
1321f40b411bSPeter Brune    Input Parameters:
132278bcb3b5SPeter Brune .  linesearch - linesearch context
1323f40b411bSPeter Brune 
1324f40b411bSPeter Brune    Output Parameters:
1325f40b411bSPeter Brune +  xnorm - The norm of the current solution
1326f40b411bSPeter Brune .  fnorm - The norm of the current function
1327f40b411bSPeter Brune -  ynorm - The norm of the current update
1328f40b411bSPeter Brune 
1329cd7522eaSPeter Brune    Notes:
1330cd7522eaSPeter Brune    This function is mainly called from SNES implementations.
1331cd7522eaSPeter Brune 
133278bcb3b5SPeter Brune    Level: developer
1333f40b411bSPeter Brune 
1334f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms() SNESLineSearchGetVecs()
1335f40b411bSPeter Brune @*/
1336f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm)
1337bf7f4e0aSPeter Brune {
1338bf7f4e0aSPeter Brune   PetscFunctionBegin;
1339f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1340f5af7f23SKarl Rupp   if (xnorm) *xnorm = linesearch->xnorm;
1341f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
1342f5af7f23SKarl Rupp   if (ynorm) *ynorm = linesearch->ynorm;
1343bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1344bf7f4e0aSPeter Brune }
1345bf7f4e0aSPeter Brune 
1346f40b411bSPeter Brune /*@
1347f1c6b773SPeter Brune    SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F.
1348f40b411bSPeter Brune 
1349f40b411bSPeter Brune    Input Parameters:
135078bcb3b5SPeter Brune +  linesearch - linesearch context
1351f40b411bSPeter Brune .  xnorm - The norm of the current solution
1352f40b411bSPeter Brune .  fnorm - The norm of the current function
1353f40b411bSPeter Brune -  ynorm - The norm of the current update
1354f40b411bSPeter Brune 
135578bcb3b5SPeter Brune    Level: advanced
1356f40b411bSPeter Brune 
1357f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1358f40b411bSPeter Brune @*/
1359f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
13606a388c36SPeter Brune {
13616a388c36SPeter Brune   PetscFunctionBegin;
1362f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
13636a388c36SPeter Brune   linesearch->xnorm = xnorm;
13646a388c36SPeter Brune   linesearch->fnorm = fnorm;
13656a388c36SPeter Brune   linesearch->ynorm = ynorm;
13666a388c36SPeter Brune   PetscFunctionReturn(0);
13676a388c36SPeter Brune }
13686a388c36SPeter Brune 
1369f40b411bSPeter Brune /*@
1370f1c6b773SPeter Brune    SNESLineSearchComputeNorms - Computes the norms of X, F, and Y.
1371f40b411bSPeter Brune 
1372f40b411bSPeter Brune    Input Parameters:
137378bcb3b5SPeter Brune .  linesearch - linesearch context
1374f40b411bSPeter Brune 
1375f40b411bSPeter Brune    Options Database Keys:
13768e557f58SPeter Brune .   -snes_linesearch_norms - turn norm computation on or off
1377f40b411bSPeter Brune 
1378f40b411bSPeter Brune    Level: intermediate
1379f40b411bSPeter Brune 
138078bcb3b5SPeter Brune .seealso: SNESLineSearchGetNorms, SNESLineSearchSetNorms(), SNESLineSearchSetComputeNorms()
1381f40b411bSPeter Brune @*/
1382f1c6b773SPeter Brune PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
13836a388c36SPeter Brune {
13846a388c36SPeter Brune   PetscErrorCode ierr;
13859bd66eb0SPeter Brune   SNES           snes;
1386bf388a1fSBarry Smith 
13876a388c36SPeter Brune   PetscFunctionBegin;
13886a388c36SPeter Brune   if (linesearch->norms) {
13899bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
1390f1c6b773SPeter Brune       ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
13919bd66eb0SPeter Brune       ierr = VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
13929bd66eb0SPeter Brune       ierr = VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
13939bd66eb0SPeter Brune       ierr = (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);CHKERRQ(ierr);
13949bd66eb0SPeter Brune     } else {
13956a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_func,   NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
13966a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_sol,    NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
13976a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
13986a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_func,     NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
13996a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_sol,      NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
14006a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_update,   NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
14016a388c36SPeter Brune     }
14029bd66eb0SPeter Brune   }
14036a388c36SPeter Brune   PetscFunctionReturn(0);
14046a388c36SPeter Brune }
14056a388c36SPeter Brune 
14066f263ca3SPeter Brune /*@
14076f263ca3SPeter Brune    SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search.
14086f263ca3SPeter Brune 
14096f263ca3SPeter Brune    Input Parameters:
141078bcb3b5SPeter Brune +  linesearch  - linesearch context
141178bcb3b5SPeter Brune -  flg  - indicates whether or not to compute norms
14126f263ca3SPeter Brune 
14136f263ca3SPeter Brune    Options Database Keys:
14141f8196a2SJed Brown .   -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic linesearch
14156f263ca3SPeter Brune 
14166f263ca3SPeter Brune    Notes:
14171f8196a2SJed Brown    This is most relevant to the SNESLINESEARCHBASIC line search type since most line searches have a stopping criteria involving the norm.
14186f263ca3SPeter Brune 
14196f263ca3SPeter Brune    Level: intermediate
14206f263ca3SPeter Brune 
14211a4f838cSPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetNorms(), SNESLineSearchComputeNorms(), SNESLINESEARCHBASIC
14226f263ca3SPeter Brune @*/
14236f263ca3SPeter Brune PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg)
14246f263ca3SPeter Brune {
14256f263ca3SPeter Brune   PetscFunctionBegin;
14266f263ca3SPeter Brune   linesearch->norms = flg;
14276f263ca3SPeter Brune   PetscFunctionReturn(0);
14286f263ca3SPeter Brune }
14296f263ca3SPeter Brune 
1430f40b411bSPeter Brune /*@
1431f1c6b773SPeter Brune    SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context
1432f40b411bSPeter Brune 
1433f40b411bSPeter Brune    Input Parameters:
143478bcb3b5SPeter Brune .  linesearch - linesearch context
1435f40b411bSPeter Brune 
1436f40b411bSPeter Brune    Output Parameters:
14376232e825SPeter Brune +  X - Solution vector
14386232e825SPeter Brune .  F - Function vector
14396232e825SPeter Brune .  Y - Search direction vector
14406232e825SPeter Brune .  W - Solution work vector
14416232e825SPeter Brune -  G - Function work vector
14426232e825SPeter Brune 
14437bba9028SPeter Brune    Notes:
14447bba9028SPeter Brune    At the beginning of a line search application, X should contain a
14456232e825SPeter Brune    solution and the vector F the function computed at X.  At the end of the
14466232e825SPeter Brune    line search application, X should contain the new solution, and F the
14476232e825SPeter Brune    function evaluated at the new solution.
1448f40b411bSPeter Brune 
14492a7a6963SBarry Smith    These vectors are owned by the SNESLineSearch and should not be destroyed by the caller
14502a7a6963SBarry Smith 
145178bcb3b5SPeter Brune    Level: advanced
1452f40b411bSPeter Brune 
1453f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1454f40b411bSPeter Brune @*/
1455bf388a1fSBarry Smith PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G)
1456bf388a1fSBarry Smith {
14576a388c36SPeter Brune   PetscFunctionBegin;
1458f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
14596a388c36SPeter Brune   if (X) {
14606a388c36SPeter Brune     PetscValidPointer(X, 2);
14616a388c36SPeter Brune     *X = linesearch->vec_sol;
14626a388c36SPeter Brune   }
14636a388c36SPeter Brune   if (F) {
14646a388c36SPeter Brune     PetscValidPointer(F, 3);
14656a388c36SPeter Brune     *F = linesearch->vec_func;
14666a388c36SPeter Brune   }
14676a388c36SPeter Brune   if (Y) {
14686a388c36SPeter Brune     PetscValidPointer(Y, 4);
14696a388c36SPeter Brune     *Y = linesearch->vec_update;
14706a388c36SPeter Brune   }
14716a388c36SPeter Brune   if (W) {
14726a388c36SPeter Brune     PetscValidPointer(W, 5);
14736a388c36SPeter Brune     *W = linesearch->vec_sol_new;
14746a388c36SPeter Brune   }
14756a388c36SPeter Brune   if (G) {
14766a388c36SPeter Brune     PetscValidPointer(G, 6);
14776a388c36SPeter Brune     *G = linesearch->vec_func_new;
14786a388c36SPeter Brune   }
14796a388c36SPeter Brune   PetscFunctionReturn(0);
14806a388c36SPeter Brune }
14816a388c36SPeter Brune 
1482f40b411bSPeter Brune /*@
1483f1c6b773SPeter Brune    SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context
1484f40b411bSPeter Brune 
1485f40b411bSPeter Brune    Input Parameters:
148678bcb3b5SPeter Brune +  linesearch - linesearch context
14876232e825SPeter Brune .  X - Solution vector
14886232e825SPeter Brune .  F - Function vector
14896232e825SPeter Brune .  Y - Search direction vector
14906232e825SPeter Brune .  W - Solution work vector
14916232e825SPeter Brune -  G - Function work vector
1492f40b411bSPeter Brune 
149378bcb3b5SPeter Brune    Level: advanced
1494f40b411bSPeter Brune 
1495f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms(), SNESLineSearchGetVecs()
1496f40b411bSPeter Brune @*/
1497bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G)
1498bf388a1fSBarry Smith {
14996a388c36SPeter Brune   PetscFunctionBegin;
1500f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
15016a388c36SPeter Brune   if (X) {
15026a388c36SPeter Brune     PetscValidHeaderSpecific(X,VEC_CLASSID,2);
15036a388c36SPeter Brune     linesearch->vec_sol = X;
15046a388c36SPeter Brune   }
15056a388c36SPeter Brune   if (F) {
15066a388c36SPeter Brune     PetscValidHeaderSpecific(F,VEC_CLASSID,3);
15076a388c36SPeter Brune     linesearch->vec_func = F;
15086a388c36SPeter Brune   }
15096a388c36SPeter Brune   if (Y) {
15106a388c36SPeter Brune     PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
15116a388c36SPeter Brune     linesearch->vec_update = Y;
15126a388c36SPeter Brune   }
15136a388c36SPeter Brune   if (W) {
15146a388c36SPeter Brune     PetscValidHeaderSpecific(W,VEC_CLASSID,5);
15156a388c36SPeter Brune     linesearch->vec_sol_new = W;
15166a388c36SPeter Brune   }
15176a388c36SPeter Brune   if (G) {
15186a388c36SPeter Brune     PetscValidHeaderSpecific(G,VEC_CLASSID,6);
15196a388c36SPeter Brune     linesearch->vec_func_new = G;
15206a388c36SPeter Brune   }
15216a388c36SPeter Brune   PetscFunctionReturn(0);
15226a388c36SPeter Brune }
15236a388c36SPeter Brune 
1524e7058c64SPeter Brune /*@C
1525f1c6b773SPeter Brune    SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1526e7058c64SPeter Brune    SNES options in the database.
1527e7058c64SPeter Brune 
1528cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
1529e7058c64SPeter Brune 
1530e7058c64SPeter Brune    Input Parameters:
1531e7058c64SPeter Brune +  snes - the SNES context
1532e7058c64SPeter Brune -  prefix - the prefix to prepend to all option names
1533e7058c64SPeter Brune 
1534e7058c64SPeter Brune    Notes:
1535e7058c64SPeter Brune    A hyphen (-) must NOT be given at the beginning of the prefix name.
1536e7058c64SPeter Brune    The first character of all runtime options is AUTOMATICALLY the hyphen.
1537e7058c64SPeter Brune 
1538e7058c64SPeter Brune    Level: advanced
1539e7058c64SPeter Brune 
1540e7058c64SPeter Brune .seealso: SNESGetOptionsPrefix()
1541e7058c64SPeter Brune @*/
1542f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[])
1543e7058c64SPeter Brune {
1544e7058c64SPeter Brune   PetscErrorCode ierr;
1545e7058c64SPeter Brune 
1546e7058c64SPeter Brune   PetscFunctionBegin;
1547f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1548e7058c64SPeter Brune   ierr = PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1549e7058c64SPeter Brune   PetscFunctionReturn(0);
1550e7058c64SPeter Brune }
1551e7058c64SPeter Brune 
1552e7058c64SPeter Brune /*@C
1553f1c6b773SPeter Brune    SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all
1554f1c6b773SPeter Brune    SNESLineSearch options in the database.
1555e7058c64SPeter Brune 
1556e7058c64SPeter Brune    Not Collective
1557e7058c64SPeter Brune 
1558e7058c64SPeter Brune    Input Parameter:
1559cd7522eaSPeter Brune .  linesearch - the SNESLineSearch context
1560e7058c64SPeter Brune 
1561e7058c64SPeter Brune    Output Parameter:
1562e7058c64SPeter Brune .  prefix - pointer to the prefix string used
1563e7058c64SPeter Brune 
15648e557f58SPeter Brune    Notes:
15658e557f58SPeter Brune    On the fortran side, the user should pass in a string 'prefix' of
1566e7058c64SPeter Brune    sufficient length to hold the prefix.
1567e7058c64SPeter Brune 
1568e7058c64SPeter Brune    Level: advanced
1569e7058c64SPeter Brune 
1570e7058c64SPeter Brune .seealso: SNESAppendOptionsPrefix()
1571e7058c64SPeter Brune @*/
1572f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[])
1573e7058c64SPeter Brune {
1574e7058c64SPeter Brune   PetscErrorCode ierr;
1575e7058c64SPeter Brune 
1576e7058c64SPeter Brune   PetscFunctionBegin;
1577f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1578e7058c64SPeter Brune   ierr = PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1579e7058c64SPeter Brune   PetscFunctionReturn(0);
1580e7058c64SPeter Brune }
1581bf7f4e0aSPeter Brune 
15828d359177SBarry Smith /*@C
1583fa0ddf94SBarry Smith    SNESLineSearchSetWorkVecs - Gets work vectors for the line search.
1584f40b411bSPeter Brune 
1585f40b411bSPeter Brune    Input Parameter:
1586f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
1587f40b411bSPeter Brune -  nwork - the number of work vectors
1588f40b411bSPeter Brune 
1589f40b411bSPeter Brune    Level: developer
1590f40b411bSPeter Brune 
1591fa0ddf94SBarry Smith .seealso: SNESSetWorkVecs()
1592f40b411bSPeter Brune @*/
1593fa0ddf94SBarry Smith PetscErrorCode  SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork)
1594bf7f4e0aSPeter Brune {
1595bf7f4e0aSPeter Brune   PetscErrorCode ierr;
1596bf388a1fSBarry Smith 
1597bf7f4e0aSPeter Brune   PetscFunctionBegin;
1598bf7f4e0aSPeter Brune   if (linesearch->vec_sol) {
1599bf7f4e0aSPeter Brune     ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr);
16008d359177SBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
1601bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1602bf7f4e0aSPeter Brune }
1603bf7f4e0aSPeter Brune 
1604f40b411bSPeter Brune /*@
1605422a814eSBarry Smith    SNESLineSearchGetReason - Gets the success/failure status of the last line search application
1606f40b411bSPeter Brune 
1607f40b411bSPeter Brune    Input Parameters:
160878bcb3b5SPeter Brune .  linesearch - linesearch context
1609f40b411bSPeter Brune 
1610f40b411bSPeter Brune    Output Parameters:
1611422a814eSBarry Smith .  result - The success or failure status
1612f40b411bSPeter Brune 
1613cd7522eaSPeter Brune    Notes:
1614c98378a5SBarry Smith    This is typically called after SNESLineSearchApply() in order to determine if the line-search failed
1615cd7522eaSPeter Brune    (and set the SNES convergence accordingly).
1616cd7522eaSPeter Brune 
1617f40b411bSPeter Brune    Level: intermediate
1618f40b411bSPeter Brune 
1619422a814eSBarry Smith .seealso: SNESLineSearchSetReason(), SNESLineSearchReason
1620f40b411bSPeter Brune @*/
1621422a814eSBarry Smith PetscErrorCode  SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result)
1622bf7f4e0aSPeter Brune {
1623bf7f4e0aSPeter Brune   PetscFunctionBegin;
1624f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1625422a814eSBarry Smith   PetscValidPointer(result, 2);
1626422a814eSBarry Smith   *result = linesearch->result;
1627bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1628bf7f4e0aSPeter Brune }
1629bf7f4e0aSPeter Brune 
1630f40b411bSPeter Brune /*@
1631422a814eSBarry Smith    SNESLineSearchSetReason - Sets the success/failure status of the last line search application
1632f40b411bSPeter Brune 
1633f40b411bSPeter Brune    Input Parameters:
163478bcb3b5SPeter Brune +  linesearch - linesearch context
1635422a814eSBarry Smith -  result - The success or failure status
1636f40b411bSPeter Brune 
1637cd7522eaSPeter Brune    Notes:
1638cd7522eaSPeter Brune    This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set
1639cd7522eaSPeter Brune    the success or failure of the line search method.
1640cd7522eaSPeter Brune 
164178bcb3b5SPeter Brune    Level: developer
1642f40b411bSPeter Brune 
1643422a814eSBarry Smith .seealso: SNESLineSearchGetSResult()
1644f40b411bSPeter Brune @*/
1645422a814eSBarry Smith PetscErrorCode  SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result)
16466a388c36SPeter Brune {
16476a388c36SPeter Brune   PetscFunctionBegin;
16485fd66863SKarl Rupp   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1649422a814eSBarry Smith   linesearch->result = result;
16506a388c36SPeter Brune   PetscFunctionReturn(0);
16516a388c36SPeter Brune }
16526a388c36SPeter Brune 
16539bd66eb0SPeter Brune /*@C
1654f1c6b773SPeter Brune    SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.
16559bd66eb0SPeter Brune 
16569bd66eb0SPeter Brune    Input Parameters:
16579bd66eb0SPeter Brune +  snes - nonlinear context obtained from SNESCreate()
16589bd66eb0SPeter Brune .  projectfunc - function for projecting the function to the bounds
16599bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
16609bd66eb0SPeter Brune 
16619bd66eb0SPeter Brune    Logically Collective on SNES
16629bd66eb0SPeter Brune 
16639bd66eb0SPeter Brune    Calling sequence of projectfunc:
16649bd66eb0SPeter Brune .vb
16659bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X)
16669bd66eb0SPeter Brune .ve
16679bd66eb0SPeter Brune 
16689bd66eb0SPeter Brune     Input parameters for projectfunc:
16699bd66eb0SPeter Brune +   snes - nonlinear context
16709bd66eb0SPeter Brune -   X - current solution
16719bd66eb0SPeter Brune 
1672cd7522eaSPeter Brune     Output parameters for projectfunc:
16739bd66eb0SPeter Brune .   X - Projected solution
16749bd66eb0SPeter Brune 
16759bd66eb0SPeter Brune    Calling sequence of normfunc:
16769bd66eb0SPeter Brune .vb
16779bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm)
16789bd66eb0SPeter Brune .ve
16799bd66eb0SPeter Brune 
1680cd7522eaSPeter Brune     Input parameters for normfunc:
16819bd66eb0SPeter Brune +   snes - nonlinear context
16829bd66eb0SPeter Brune .   X - current solution
16839bd66eb0SPeter Brune -   F - current residual
16849bd66eb0SPeter Brune 
1685cd7522eaSPeter Brune     Output parameters for normfunc:
16869bd66eb0SPeter Brune .   fnorm - VI-specific norm of the function
16879bd66eb0SPeter Brune 
1688cd7522eaSPeter Brune     Notes:
1689cd7522eaSPeter Brune     The VI solvers require projection of the solution to the feasible set.  projectfunc should implement this.
1690cd7522eaSPeter Brune 
1691cd7522eaSPeter Brune     The VI solvers require special evaluation of the function norm such that the norm is only calculated
1692cd7522eaSPeter Brune     on the inactive set.  This should be implemented by normfunc.
16939bd66eb0SPeter Brune 
16949bd66eb0SPeter Brune     Level: developer
16959bd66eb0SPeter Brune 
1696f1c6b773SPeter Brune .seealso: SNESLineSearchGetVIFunctions(), SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck()
16979bd66eb0SPeter Brune @*/
169825acbd8eSLisandro Dalcin PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc)
16999bd66eb0SPeter Brune {
17009bd66eb0SPeter Brune   PetscFunctionBegin;
1701f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
17029bd66eb0SPeter Brune   if (projectfunc) linesearch->ops->viproject = projectfunc;
17039bd66eb0SPeter Brune   if (normfunc) linesearch->ops->vinorm = normfunc;
17049bd66eb0SPeter Brune   PetscFunctionReturn(0);
17059bd66eb0SPeter Brune }
17069bd66eb0SPeter Brune 
17079bd66eb0SPeter Brune /*@C
1708f1c6b773SPeter Brune    SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.
17099bd66eb0SPeter Brune 
17109bd66eb0SPeter Brune    Input Parameters:
1711907376e6SBarry Smith .  linesearch - the line search context, obtain with SNESGetLineSearch()
17129bd66eb0SPeter Brune 
17139bd66eb0SPeter Brune    Output Parameters:
17149bd66eb0SPeter Brune +  projectfunc - function for projecting the function to the bounds
17159bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
17169bd66eb0SPeter Brune 
17179bd66eb0SPeter Brune    Logically Collective on SNES
17189bd66eb0SPeter Brune 
17199bd66eb0SPeter Brune     Level: developer
17209bd66eb0SPeter Brune 
1721f1c6b773SPeter Brune .seealso: SNESLineSearchSetVIFunctions(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
17229bd66eb0SPeter Brune @*/
172325acbd8eSLisandro Dalcin PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc)
17249bd66eb0SPeter Brune {
17259bd66eb0SPeter Brune   PetscFunctionBegin;
17269bd66eb0SPeter Brune   if (projectfunc) *projectfunc = linesearch->ops->viproject;
17279bd66eb0SPeter Brune   if (normfunc) *normfunc = linesearch->ops->vinorm;
17289bd66eb0SPeter Brune   PetscFunctionReturn(0);
17299bd66eb0SPeter Brune }
17309bd66eb0SPeter Brune 
1731bf7f4e0aSPeter Brune /*@C
17321c84c290SBarry Smith   SNESLineSearchRegister - See SNESLineSearchRegister()
1733bf7f4e0aSPeter Brune 
1734bf7f4e0aSPeter Brune   Level: advanced
1735bf7f4e0aSPeter Brune @*/
1736bdf89e91SBarry Smith PetscErrorCode  SNESLineSearchRegister(const char sname[],PetscErrorCode (*function)(SNESLineSearch))
1737bf7f4e0aSPeter Brune {
1738bf7f4e0aSPeter Brune   PetscErrorCode ierr;
1739bf7f4e0aSPeter Brune 
1740bf7f4e0aSPeter Brune   PetscFunctionBegin;
17411d36bdfdSBarry Smith   ierr = SNESInitializePackage();CHKERRQ(ierr);
1742a240a19fSJed Brown   ierr = PetscFunctionListAdd(&SNESLineSearchList,sname,function);CHKERRQ(ierr);
1743bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1744bf7f4e0aSPeter Brune }
1745