xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision 534a8f05a7a8aff70dd8cfd53d9cd834400a8dbf)
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;
209bf7f4e0aSPeter Brune   *outlinesearch           = linesearch;
210bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
211bf7f4e0aSPeter Brune }
212bf7f4e0aSPeter Brune 
213f40b411bSPeter Brune /*@
21478bcb3b5SPeter Brune    SNESLineSearchSetUp - Prepares the line search for being applied by allocating
21578bcb3b5SPeter Brune    any required vectors.
216f40b411bSPeter Brune 
217cd7522eaSPeter Brune    Collective on SNESLineSearch
218f40b411bSPeter Brune 
219f40b411bSPeter Brune    Input Parameters:
220f40b411bSPeter Brune .  linesearch - The LineSearch instance.
221f40b411bSPeter Brune 
222cd7522eaSPeter Brune    Notes:
223f190f2fcSBarry Smith    For most cases, this needn't be called by users or outside of SNESLineSearchApply().
224cd7522eaSPeter Brune    The only current case where this is called outside of this is for the VI
22578bcb3b5SPeter Brune    solvers, which modify the solution and work vectors before the first call
226cd7522eaSPeter Brune    of SNESLineSearchApply, requiring the SNESLineSearch work vectors to be
227cd7522eaSPeter Brune    allocated upfront.
228cd7522eaSPeter Brune 
22978bcb3b5SPeter Brune    Level: advanced
230f40b411bSPeter Brune 
23184238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchReset()
232f40b411bSPeter Brune @*/
233f40b411bSPeter Brune 
234bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch)
235bf388a1fSBarry Smith {
236bf7f4e0aSPeter Brune   PetscErrorCode ierr;
237bf388a1fSBarry Smith 
238bf7f4e0aSPeter Brune   PetscFunctionBegin;
239bf7f4e0aSPeter Brune   if (!((PetscObject)linesearch)->type_name) {
2401a4f838cSPeter Brune     ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr);
241bf7f4e0aSPeter Brune   }
242bf7f4e0aSPeter Brune   if (!linesearch->setupcalled) {
2439bd66eb0SPeter Brune     if (!linesearch->vec_sol_new) {
244bf7f4e0aSPeter Brune       ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new);CHKERRQ(ierr);
2459bd66eb0SPeter Brune     }
2469bd66eb0SPeter Brune     if (!linesearch->vec_func_new) {
2479bd66eb0SPeter Brune       ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new);CHKERRQ(ierr);
2489bd66eb0SPeter Brune     }
249bf7f4e0aSPeter Brune     if (linesearch->ops->setup) {
250bf7f4e0aSPeter Brune       ierr = (*linesearch->ops->setup)(linesearch);CHKERRQ(ierr);
251bf7f4e0aSPeter Brune     }
252ed07d7d7SPeter Brune     if (!linesearch->ops->snesfunc) {ierr = SNESLineSearchSetFunction(linesearch,SNESComputeFunction);CHKERRQ(ierr);}
253bf7f4e0aSPeter Brune     linesearch->lambda      = linesearch->damping;
254bf7f4e0aSPeter Brune     linesearch->setupcalled = PETSC_TRUE;
255bf7f4e0aSPeter Brune   }
256bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
257bf7f4e0aSPeter Brune }
258bf7f4e0aSPeter Brune 
259f40b411bSPeter Brune 
260f40b411bSPeter Brune /*@
261f190f2fcSBarry Smith    SNESLineSearchReset - Undoes the SNESLineSearchSetUp() and deletes any Vecs or Mats allocated by the line search.
262f40b411bSPeter Brune 
263f1c6b773SPeter Brune    Collective on SNESLineSearch
264f40b411bSPeter Brune 
265f40b411bSPeter Brune    Input Parameters:
266f40b411bSPeter Brune .  linesearch - The LineSearch instance.
267f40b411bSPeter Brune 
26895452b02SPatrick Sanan    Notes:
26995452b02SPatrick Sanan     Usually only called by SNESReset()
270f190f2fcSBarry Smith 
271f190f2fcSBarry Smith    Level: developer
272f40b411bSPeter Brune 
27384238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchSetUp()
274f40b411bSPeter Brune @*/
275f40b411bSPeter Brune 
276bf388a1fSBarry Smith PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch)
277bf388a1fSBarry Smith {
278bf7f4e0aSPeter Brune   PetscErrorCode ierr;
279bf388a1fSBarry Smith 
280bf7f4e0aSPeter Brune   PetscFunctionBegin;
281f5af7f23SKarl Rupp   if (linesearch->ops->reset) (*linesearch->ops->reset)(linesearch);
282f5af7f23SKarl Rupp 
283bf7f4e0aSPeter Brune   ierr = VecDestroy(&linesearch->vec_sol_new);CHKERRQ(ierr);
284bf7f4e0aSPeter Brune   ierr = VecDestroy(&linesearch->vec_func_new);CHKERRQ(ierr);
285bf7f4e0aSPeter Brune 
286bf7f4e0aSPeter Brune   ierr = VecDestroyVecs(linesearch->nwork, &linesearch->work);CHKERRQ(ierr);
287f5af7f23SKarl Rupp 
288bf7f4e0aSPeter Brune   linesearch->nwork       = 0;
289bf7f4e0aSPeter Brune   linesearch->setupcalled = PETSC_FALSE;
290bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
291bf7f4e0aSPeter Brune }
292bf7f4e0aSPeter Brune 
293ed07d7d7SPeter Brune /*@C
294f190f2fcSBarry Smith    SNESLineSearchSetFunction - Sets the function evaluation used by the SNES line search
295ed07d7d7SPeter Brune 
296ed07d7d7SPeter Brune    Input Parameters:
297ed07d7d7SPeter Brune .  linesearch - the SNESLineSearch context
298f190f2fcSBarry Smith +  func       - function evaluation routine
299ed07d7d7SPeter Brune 
300ed07d7d7SPeter Brune    Level: developer
301ed07d7d7SPeter Brune 
30295452b02SPatrick Sanan    Notes:
30395452b02SPatrick Sanan     This is used internally by PETSc and not called by users
304f190f2fcSBarry Smith 
30584238204SBarry Smith .seealso: SNESGetLineSearch(), SNESSetFunction()
306ed07d7d7SPeter Brune @*/
307ed07d7d7SPeter Brune PetscErrorCode  SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES,Vec,Vec))
308ed07d7d7SPeter Brune {
309ed07d7d7SPeter Brune   PetscFunctionBegin;
310ed07d7d7SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
311ed07d7d7SPeter Brune   linesearch->ops->snesfunc = func;
312ed07d7d7SPeter Brune   PetscFunctionReturn(0);
313ed07d7d7SPeter Brune }
314ed07d7d7SPeter Brune 
31586d74e61SPeter Brune /*@C
316f190f2fcSBarry Smith    SNESLineSearchSetPreCheck - Sets a user function that is called after the initial search direction has been computed but
317df3898eeSBarry Smith          before the line search routine has been applied. Allows the user to adjust the result of (usually a linear solve) that
318f190f2fcSBarry Smith          determined the search direction.
31986d74e61SPeter Brune 
320f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
32186d74e61SPeter Brune 
32286d74e61SPeter Brune    Input Parameters:
323f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
32484238204SBarry Smith .  func - [optional] function evaluation routine, see SNESLineSearchPreCheck() for the calling sequence
325f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
32686d74e61SPeter Brune 
32786d74e61SPeter Brune    Level: intermediate
32886d74e61SPeter Brune 
32984238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
33086d74e61SPeter Brune @*/
331f190f2fcSBarry Smith PetscErrorCode  SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void *ctx)
33286d74e61SPeter Brune {
3339bd66eb0SPeter Brune   PetscFunctionBegin;
334f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
335f190f2fcSBarry Smith   if (func) linesearch->ops->precheck = func;
33686d74e61SPeter Brune   if (ctx) linesearch->precheckctx = ctx;
33786d74e61SPeter Brune   PetscFunctionReturn(0);
33886d74e61SPeter Brune }
33986d74e61SPeter Brune 
34086d74e61SPeter Brune /*@C
34153deb899SBarry Smith    SNESLineSearchGetPreCheck - Gets the pre-check function for the line search routine.
34286d74e61SPeter Brune 
34386d74e61SPeter Brune    Input Parameters:
344f1c6b773SPeter Brune .  linesearch - the SNESLineSearch context
34586d74e61SPeter Brune 
34686d74e61SPeter Brune    Output Parameters:
34784238204SBarry Smith +  func       - [optional] function evaluation routine, see SNESLineSearchPreCheck() for calling sequence
348f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
34986d74e61SPeter Brune 
35086d74e61SPeter Brune    Level: intermediate
35186d74e61SPeter Brune 
35284238204SBarry Smith .seealso: SNESGetLineSearch(), SNESGetLineSearch(), SNESLineSearchPreCheck(), SNESLineSearchGetPostCheck(), SNESLineSearchSetPreCheck(), SNESLineSearchSetPostCheck()
35386d74e61SPeter Brune @*/
354f190f2fcSBarry Smith PetscErrorCode  SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void **ctx)
35586d74e61SPeter Brune {
3569bd66eb0SPeter Brune   PetscFunctionBegin;
357f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
358f190f2fcSBarry Smith   if (func) *func = linesearch->ops->precheck;
35986d74e61SPeter Brune   if (ctx) *ctx = linesearch->precheckctx;
36086d74e61SPeter Brune   PetscFunctionReturn(0);
36186d74e61SPeter Brune }
36286d74e61SPeter Brune 
36386d74e61SPeter Brune /*@C
364f190f2fcSBarry Smith    SNESLineSearchSetPostCheck - Sets a user function that is called after the line search has been applied to determine the step
365f190f2fcSBarry Smith        direction and length. Allows the user a chance to change or override the decision of the line search routine
36686d74e61SPeter Brune 
367f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
36886d74e61SPeter Brune 
36986d74e61SPeter Brune    Input Parameters:
370f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
37184238204SBarry Smith .  func - [optional] function evaluation routine, see SNESLineSearchPostCheck()  for the calling sequence
372f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
37386d74e61SPeter Brune 
37486d74e61SPeter Brune    Level: intermediate
37586d74e61SPeter Brune 
37684238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchPostCheck(), SNESLineSearchSetPreCheck(), SNESLineSearchGetPreCheck(), SNESLineSearchGetPostCheck()
37786d74e61SPeter Brune @*/
378f190f2fcSBarry Smith PetscErrorCode  SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void *ctx)
37986d74e61SPeter Brune {
38086d74e61SPeter Brune   PetscFunctionBegin;
381f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
382f190f2fcSBarry Smith   if (func) linesearch->ops->postcheck = func;
38386d74e61SPeter Brune   if (ctx) linesearch->postcheckctx = ctx;
38486d74e61SPeter Brune   PetscFunctionReturn(0);
38586d74e61SPeter Brune }
38686d74e61SPeter Brune 
38786d74e61SPeter Brune /*@C
388f1c6b773SPeter Brune    SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine.
38986d74e61SPeter Brune 
39086d74e61SPeter Brune    Input Parameters:
391f1c6b773SPeter Brune .  linesearch - the SNESLineSearch context
39286d74e61SPeter Brune 
39386d74e61SPeter Brune    Output Parameters:
39484238204SBarry Smith +  func - [optional] function evaluation routine, see for the calling sequence SNESLineSearchPostCheck()
395f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
39686d74e61SPeter Brune 
39786d74e61SPeter Brune    Level: intermediate
39886d74e61SPeter Brune 
39984238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchPostCheck(), SNESLineSearchSetPreCheck()
40086d74e61SPeter Brune @*/
401f190f2fcSBarry Smith PetscErrorCode  SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void **ctx)
40286d74e61SPeter Brune {
4039bd66eb0SPeter Brune   PetscFunctionBegin;
404f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
405f190f2fcSBarry Smith   if (func) *func = linesearch->ops->postcheck;
40686d74e61SPeter Brune   if (ctx) *ctx = linesearch->postcheckctx;
40786d74e61SPeter Brune   PetscFunctionReturn(0);
40886d74e61SPeter Brune }
40986d74e61SPeter Brune 
410f40b411bSPeter Brune /*@
411f1c6b773SPeter Brune    SNESLineSearchPreCheck - Prepares the line search for being applied.
412f40b411bSPeter Brune 
413cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
414f40b411bSPeter Brune 
415f40b411bSPeter Brune    Input Parameters:
4167b1df9c1SPeter Brune +  linesearch - The linesearch instance.
4177b1df9c1SPeter Brune .  X - The current solution
4187b1df9c1SPeter Brune -  Y - The step direction
419f40b411bSPeter Brune 
420f40b411bSPeter Brune    Output Parameters:
4218e557f58SPeter Brune .  changed - Indicator that the precheck routine has changed anything
422f40b411bSPeter Brune 
423f190f2fcSBarry Smith    Level: developer
424f40b411bSPeter Brune 
42584238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchPostCheck(), SNESLineSearchSetPreCheck(), SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck()
426f40b411bSPeter Brune @*/
4277b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed)
428bf7f4e0aSPeter Brune {
429bf7f4e0aSPeter Brune   PetscErrorCode ierr;
4305fd66863SKarl Rupp 
431bf7f4e0aSPeter Brune   PetscFunctionBegin;
432bf7f4e0aSPeter Brune   *changed = PETSC_FALSE;
4336b2b7091SBarry Smith   if (linesearch->ops->precheck) {
4346b2b7091SBarry Smith     ierr = (*linesearch->ops->precheck)(linesearch, X, Y, changed, linesearch->precheckctx);CHKERRQ(ierr);
43538bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed,4);
436bf7f4e0aSPeter Brune   }
437bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
438bf7f4e0aSPeter Brune }
439bf7f4e0aSPeter Brune 
440f40b411bSPeter Brune /*@
441f1c6b773SPeter Brune    SNESLineSearchPostCheck - Prepares the line search for being applied.
442f40b411bSPeter Brune 
443cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
444f40b411bSPeter Brune 
445f40b411bSPeter Brune    Input Parameters:
4467b1df9c1SPeter Brune +  linesearch - The linesearch context
4477b1df9c1SPeter Brune .  X - The last solution
4487b1df9c1SPeter Brune .  Y - The step direction
4497b1df9c1SPeter Brune -  W - The updated solution, W = X + lambda*Y for some lambda
450f40b411bSPeter Brune 
451f40b411bSPeter Brune    Output Parameters:
45278bcb3b5SPeter Brune +  changed_Y - Indicator if the direction Y has been changed.
45378bcb3b5SPeter Brune -  changed_W - Indicator if the new candidate solution W has been changed.
454f40b411bSPeter Brune 
455f190f2fcSBarry Smith    Level: developer
456f40b411bSPeter Brune 
45784238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck(), SNESLineSearchSetPrecheck(), SNESLineSearchGetPrecheck()
458f40b411bSPeter Brune @*/
4597b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch,Vec X,Vec Y,Vec W,PetscBool *changed_Y,PetscBool *changed_W)
460bf7f4e0aSPeter Brune {
461bf7f4e0aSPeter Brune   PetscErrorCode ierr;
462bf388a1fSBarry Smith 
463bf7f4e0aSPeter Brune   PetscFunctionBegin;
464bf7f4e0aSPeter Brune   *changed_Y = PETSC_FALSE;
465bf7f4e0aSPeter Brune   *changed_W = PETSC_FALSE;
4666b2b7091SBarry Smith   if (linesearch->ops->postcheck) {
4676b2b7091SBarry Smith     ierr = (*linesearch->ops->postcheck)(linesearch,X,Y,W,changed_Y,changed_W,linesearch->postcheckctx);CHKERRQ(ierr);
46838bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed_Y,5);
46938bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed_W,6);
47086d74e61SPeter Brune   }
47186d74e61SPeter Brune   PetscFunctionReturn(0);
47286d74e61SPeter Brune }
47386d74e61SPeter Brune 
47486d74e61SPeter Brune /*@C
47586d74e61SPeter Brune    SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration
47686d74e61SPeter Brune 
477cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
47886d74e61SPeter Brune 
47986d74e61SPeter Brune    Input Arguments:
48086d74e61SPeter Brune +  linesearch - linesearch context
48186d74e61SPeter Brune .  X - base state for this step
48286d74e61SPeter Brune .  Y - initial correction
483907376e6SBarry Smith -  ctx - context for this function
48486d74e61SPeter Brune 
48586d74e61SPeter Brune    Output Arguments:
48686d74e61SPeter Brune +  Y - correction, possibly modified
48786d74e61SPeter Brune -  changed - flag indicating that Y was modified
48886d74e61SPeter Brune 
48986d74e61SPeter Brune    Options Database Key:
490cd7522eaSPeter Brune +  -snes_linesearch_precheck_picard - activate this routine
491cd7522eaSPeter Brune -  -snes_linesearch_precheck_picard_angle - angle
49286d74e61SPeter Brune 
49386d74e61SPeter Brune    Level: advanced
49486d74e61SPeter Brune 
49586d74e61SPeter Brune    Notes:
49686d74e61SPeter Brune    This function should be passed to SNESLineSearchSetPreCheck()
49786d74e61SPeter Brune 
49886d74e61SPeter Brune    The justification for this method involves the linear convergence of a Picard iteration
49986d74e61SPeter Brune    so the Picard linearization should be provided in place of the "Jacobian". This correction
50086d74e61SPeter Brune    is generally not useful when using a Newton linearization.
50186d74e61SPeter Brune 
50286d74e61SPeter Brune    Reference:
50386d74e61SPeter Brune    Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology.
50486d74e61SPeter Brune 
50584238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchSetPreCheck()
50686d74e61SPeter Brune @*/
507f1c6b773SPeter Brune PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed,void *ctx)
50886d74e61SPeter Brune {
50986d74e61SPeter Brune   PetscErrorCode ierr;
51086d74e61SPeter Brune   PetscReal      angle = *(PetscReal*)linesearch->precheckctx;
51186d74e61SPeter Brune   Vec            Ylast;
51286d74e61SPeter Brune   PetscScalar    dot;
51386d74e61SPeter Brune   PetscInt       iter;
51486d74e61SPeter Brune   PetscReal      ynorm,ylastnorm,theta,angle_radians;
51586d74e61SPeter Brune   SNES           snes;
51686d74e61SPeter Brune 
51786d74e61SPeter Brune   PetscFunctionBegin;
518f1c6b773SPeter Brune   ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
51986d74e61SPeter Brune   ierr = PetscObjectQuery((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject*)&Ylast);CHKERRQ(ierr);
52086d74e61SPeter Brune   if (!Ylast) {
52186d74e61SPeter Brune     ierr = VecDuplicate(Y,&Ylast);CHKERRQ(ierr);
52286d74e61SPeter Brune     ierr = PetscObjectCompose((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject)Ylast);CHKERRQ(ierr);
52386d74e61SPeter Brune     ierr = PetscObjectDereference((PetscObject)Ylast);CHKERRQ(ierr);
52486d74e61SPeter Brune   }
52586d74e61SPeter Brune   ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr);
52686d74e61SPeter Brune   if (iter < 2) {
52786d74e61SPeter Brune     ierr     = VecCopy(Y,Ylast);CHKERRQ(ierr);
52886d74e61SPeter Brune     *changed = PETSC_FALSE;
52986d74e61SPeter Brune     PetscFunctionReturn(0);
53086d74e61SPeter Brune   }
53186d74e61SPeter Brune 
53286d74e61SPeter Brune   ierr = VecDot(Y,Ylast,&dot);CHKERRQ(ierr);
53386d74e61SPeter Brune   ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr);
53486d74e61SPeter Brune   ierr = VecNorm(Ylast,NORM_2,&ylastnorm);CHKERRQ(ierr);
53586d74e61SPeter Brune   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
536255453a1SBarry Smith   theta         = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0));
53786d74e61SPeter Brune   angle_radians = angle * PETSC_PI / 180.;
53886d74e61SPeter Brune   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
53986d74e61SPeter Brune     /* Modify the step Y */
54086d74e61SPeter Brune     PetscReal alpha,ydiffnorm;
54186d74e61SPeter Brune     ierr  = VecAXPY(Ylast,-1.0,Y);CHKERRQ(ierr);
54286d74e61SPeter Brune     ierr  = VecNorm(Ylast,NORM_2,&ydiffnorm);CHKERRQ(ierr);
543f85e2ce2SBarry Smith     alpha = (ydiffnorm > .001*ylastnorm) ? ylastnorm / ydiffnorm : 1000.0;
54486d74e61SPeter Brune     ierr  = VecCopy(Y,Ylast);CHKERRQ(ierr);
54586d74e61SPeter Brune     ierr  = VecScale(Y,alpha);CHKERRQ(ierr);
546c69d1a72SBarry 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);
547a47ec85fSBarry Smith     *changed = PETSC_TRUE;
54886d74e61SPeter Brune   } else {
549c69d1a72SBarry 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);
55086d74e61SPeter Brune     ierr     = VecCopy(Y,Ylast);CHKERRQ(ierr);
55186d74e61SPeter Brune     *changed = PETSC_FALSE;
552bf7f4e0aSPeter Brune   }
553bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
554bf7f4e0aSPeter Brune }
555bf7f4e0aSPeter Brune 
556f40b411bSPeter Brune /*@
557cd7522eaSPeter Brune    SNESLineSearchApply - Computes the line-search update.
558f40b411bSPeter Brune 
559f1c6b773SPeter Brune    Collective on SNESLineSearch
560f40b411bSPeter Brune 
561f40b411bSPeter Brune    Input Parameters:
5628e557f58SPeter Brune +  linesearch - The linesearch context
5638e557f58SPeter Brune .  X - The current solution
5648e557f58SPeter Brune .  F - The current function
5658e557f58SPeter Brune .  fnorm - The current norm
5668e557f58SPeter Brune -  Y - The search direction
567f40b411bSPeter Brune 
568f40b411bSPeter Brune    Output Parameters:
5698e557f58SPeter Brune +  X - The new solution
5708e557f58SPeter Brune .  F - The new function
5718e557f58SPeter Brune -  fnorm - The new function norm
572f40b411bSPeter Brune 
573cd7522eaSPeter Brune    Options Database Keys:
574d4c6564cSPatrick Farrell + -snes_linesearch_type - basic, bt, l2, cp, nleqerr, shell
575dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches
5761fe24845SBarry Smith . -snes_linesearch_damping - The linesearch damping parameter, default is 1.0 (no damping)
5771fe24845SBarry Smith . -snes_linesearch_norms   - Turn on/off the linesearch norms computation (SNESLineSearchSetComputeNorms())
5783c7d6663SPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess
5793c7d6663SPeter Brune - -snes_linesearch_max_it - The number of iterations for iterative line searches
580cd7522eaSPeter Brune 
581cd7522eaSPeter Brune    Notes:
582cd7522eaSPeter Brune    This is typically called from within a SNESSolve() implementation in order to
583cd7522eaSPeter Brune    help with convergence of the nonlinear method.  Various SNES types use line searches
584cd7522eaSPeter Brune    in different ways, but the overarching theme is that a line search is used to determine
585cd7522eaSPeter Brune    an optimal damping parameter of a step at each iteration of the method.  Each
58684238204SBarry Smith    application of the line search may invoke SNESComputeFunction() several times, and
587cd7522eaSPeter Brune    therefore may be fairly expensive.
588cd7522eaSPeter Brune 
589f40b411bSPeter Brune    Level: Intermediate
590f40b411bSPeter Brune 
59184238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchCreate(), SNESLineSearchPreCheck(), SNESLineSearchPostCheck(), SNESSolve(), SNESComputeFunction(), SNESLineSearchSetComputeNorms(),
5921fe24845SBarry Smith           SNESLineSearchType, SNESLineSearchSetType()
593f40b411bSPeter Brune @*/
594bf388a1fSBarry Smith PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y)
595bf388a1fSBarry Smith {
596bf7f4e0aSPeter Brune   PetscErrorCode ierr;
597bf7f4e0aSPeter Brune 
598bf388a1fSBarry Smith   PetscFunctionBegin;
599f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
600bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
601bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
602bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
603bf7f4e0aSPeter Brune 
604422a814eSBarry Smith   linesearch->result = SNES_LINESEARCH_SUCCEEDED;
605bf7f4e0aSPeter Brune 
606bf7f4e0aSPeter Brune   linesearch->vec_sol    = X;
607bf7f4e0aSPeter Brune   linesearch->vec_update = Y;
608bf7f4e0aSPeter Brune   linesearch->vec_func   = F;
609bf7f4e0aSPeter Brune 
610f1c6b773SPeter Brune   ierr = SNESLineSearchSetUp(linesearch);CHKERRQ(ierr);
611bf7f4e0aSPeter Brune 
612f5af7f23SKarl Rupp   if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */
613bf7f4e0aSPeter Brune 
614f5af7f23SKarl Rupp   if (fnorm) linesearch->fnorm = *fnorm;
615f5af7f23SKarl Rupp   else {
616bf7f4e0aSPeter Brune     ierr = VecNorm(F, NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
617bf7f4e0aSPeter Brune   }
618bf7f4e0aSPeter Brune 
61957a83faaSBarry Smith   ierr = PetscLogEventBegin(SNESLINESEARCH_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
620bf7f4e0aSPeter Brune 
621bf7f4e0aSPeter Brune   ierr = (*linesearch->ops->apply)(linesearch);CHKERRQ(ierr);
622bf7f4e0aSPeter Brune 
62357a83faaSBarry Smith   ierr = PetscLogEventEnd(SNESLINESEARCH_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
624bf7f4e0aSPeter Brune 
625f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
626bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
627bf7f4e0aSPeter Brune }
628bf7f4e0aSPeter Brune 
629f40b411bSPeter Brune /*@
630f1c6b773SPeter Brune    SNESLineSearchDestroy - Destroys the line search instance.
631f40b411bSPeter Brune 
632f1c6b773SPeter Brune    Collective on SNESLineSearch
633f40b411bSPeter Brune 
634f40b411bSPeter Brune    Input Parameters:
6358e557f58SPeter Brune .  linesearch - The linesearch context
636f40b411bSPeter Brune 
63784238204SBarry Smith    Level: developer
638f40b411bSPeter Brune 
63984238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchCreate(), SNESLineSearchReset(), SNESDestroy()
640f40b411bSPeter Brune @*/
641bf388a1fSBarry Smith PetscErrorCode SNESLineSearchDestroy(SNESLineSearch * linesearch)
642bf388a1fSBarry Smith {
643bf7f4e0aSPeter Brune   PetscErrorCode ierr;
644bf388a1fSBarry Smith 
645bf7f4e0aSPeter Brune   PetscFunctionBegin;
646bf7f4e0aSPeter Brune   if (!*linesearch) PetscFunctionReturn(0);
647f1c6b773SPeter Brune   PetscValidHeaderSpecific((*linesearch),SNESLINESEARCH_CLASSID,1);
648bf7f4e0aSPeter Brune   if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = 0; PetscFunctionReturn(0);}
649e04113cfSBarry Smith   ierr = PetscObjectSAWsViewOff((PetscObject)*linesearch);CHKERRQ(ierr);
65022d28d08SBarry Smith   ierr = SNESLineSearchReset(*linesearch);CHKERRQ(ierr);
651f5af7f23SKarl Rupp   if ((*linesearch)->ops->destroy) (*linesearch)->ops->destroy(*linesearch);
652bf7f4e0aSPeter Brune   ierr = PetscViewerDestroy(&(*linesearch)->monitor);CHKERRQ(ierr);
653dcf2fd19SBarry Smith   ierr = SNESLineSearchMonitorCancel((*linesearch));CHKERRQ(ierr);
654e7058c64SPeter Brune   ierr = PetscHeaderDestroy(linesearch);CHKERRQ(ierr);
655bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
656bf7f4e0aSPeter Brune }
657bf7f4e0aSPeter Brune 
658f40b411bSPeter Brune /*@
659dcf2fd19SBarry Smith    SNESLineSearchSetDefaultMonitor - Turns on/off printing useful information and debugging output about the line search.
660bf7f4e0aSPeter Brune 
661bf7f4e0aSPeter Brune    Input Parameters:
662dcf2fd19SBarry Smith +  linesearch - the linesearch object
663dcf2fd19SBarry Smith -  viewer - an ASCII PetscViewer or NULL to turn off monitor
664bf7f4e0aSPeter Brune 
665dcf2fd19SBarry Smith    Logically Collective on SNESLineSearch
666bf7f4e0aSPeter Brune 
667bf7f4e0aSPeter Brune    Options Database:
668dcf2fd19SBarry Smith .   -snes_linesearch_monitor [:filename] - enables the monitor
669bf7f4e0aSPeter Brune 
670bf7f4e0aSPeter Brune    Level: intermediate
671bf7f4e0aSPeter Brune 
672dcf2fd19SBarry Smith    Developer Note: This monitor is implemented differently than the other SNESLineSearchMonitors that are set with
673d12e167eSBarry Smith      SNESLineSearchMonitorSet() since it is called in many locations of the line search routines to display aspects of the
674d12e167eSBarry Smith      line search that are not visible to the other monitors.
675dcf2fd19SBarry Smith 
67684238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchGetDefaultMonitor(), PetscViewer, SNESLineSearchSetMonitor()
677bf7f4e0aSPeter Brune @*/
678dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchSetDefaultMonitor(SNESLineSearch linesearch, PetscViewer viewer)
679bf7f4e0aSPeter Brune {
680bf7f4e0aSPeter Brune   PetscErrorCode ierr;
681bf388a1fSBarry Smith 
682bf7f4e0aSPeter Brune   PetscFunctionBegin;
683dcf2fd19SBarry Smith   if (viewer) {ierr = PetscObjectReference((PetscObject)viewer);CHKERRQ(ierr);}
684bf7f4e0aSPeter Brune   ierr = PetscViewerDestroy(&linesearch->monitor);CHKERRQ(ierr);
685dcf2fd19SBarry Smith   linesearch->monitor = viewer;
686bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
687bf7f4e0aSPeter Brune }
688bf7f4e0aSPeter Brune 
689f40b411bSPeter Brune /*@
690dcf2fd19SBarry Smith    SNESLineSearchGetDefaultMonitor - Gets the PetscViewer instance for the line search monitor.
6916a388c36SPeter Brune 
692f190f2fcSBarry Smith    Input Parameter:
6938e557f58SPeter Brune .  linesearch - linesearch context
694f40b411bSPeter Brune 
695f190f2fcSBarry Smith    Output Parameter:
6968e557f58SPeter Brune .  monitor - monitor context
697f40b411bSPeter Brune 
698f40b411bSPeter Brune    Logically Collective on SNES
699f40b411bSPeter Brune 
700f40b411bSPeter Brune    Options Database Keys:
7018e557f58SPeter Brune .   -snes_linesearch_monitor - enables the monitor
702f40b411bSPeter Brune 
703f40b411bSPeter Brune    Level: intermediate
704f40b411bSPeter Brune 
70584238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchSetDefaultMonitor(), PetscViewer
706f40b411bSPeter Brune @*/
707dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchGetDefaultMonitor(SNESLineSearch linesearch, PetscViewer *monitor)
7086a388c36SPeter Brune {
7096a388c36SPeter Brune   PetscFunctionBegin;
710f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
7116a388c36SPeter Brune   if (monitor) {
7126a388c36SPeter Brune     PetscValidPointer(monitor, 2);
7136a388c36SPeter Brune     *monitor = linesearch->monitor;
7146a388c36SPeter Brune   }
7156a388c36SPeter Brune   PetscFunctionReturn(0);
7166a388c36SPeter Brune }
7176a388c36SPeter Brune 
718dcf2fd19SBarry Smith /*@C
719dcf2fd19SBarry Smith    SNESLineSearchMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
720dcf2fd19SBarry Smith 
721dcf2fd19SBarry Smith    Collective on SNESLineSearch
722dcf2fd19SBarry Smith 
723dcf2fd19SBarry Smith    Input Parameters:
724dcf2fd19SBarry Smith +  ls - LineSearch object you wish to monitor
725dcf2fd19SBarry Smith .  name - the monitor type one is seeking
726dcf2fd19SBarry Smith .  help - message indicating what monitoring is done
727dcf2fd19SBarry Smith .  manual - manual page for the monitor
728dcf2fd19SBarry Smith .  monitor - the monitor function
729dcf2fd19SBarry 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
730dcf2fd19SBarry Smith 
731dcf2fd19SBarry Smith    Level: developer
732dcf2fd19SBarry Smith 
733dcf2fd19SBarry Smith .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
734dcf2fd19SBarry Smith           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
735dcf2fd19SBarry Smith           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
736dcf2fd19SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
737dcf2fd19SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
738dcf2fd19SBarry Smith           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
739dcf2fd19SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
740dcf2fd19SBarry Smith @*/
741d12e167eSBarry Smith PetscErrorCode  SNESLineSearchMonitorSetFromOptions(SNESLineSearch ls,const char name[],const char help[], const char manual[],PetscErrorCode (*monitor)(SNESLineSearch,PetscViewerAndFormat*),PetscErrorCode (*monitorsetup)(SNESLineSearch,PetscViewerAndFormat*))
742dcf2fd19SBarry Smith {
743dcf2fd19SBarry Smith   PetscErrorCode    ierr;
744dcf2fd19SBarry Smith   PetscViewer       viewer;
745dcf2fd19SBarry Smith   PetscViewerFormat format;
746dcf2fd19SBarry Smith   PetscBool         flg;
747dcf2fd19SBarry Smith 
748dcf2fd19SBarry Smith   PetscFunctionBegin;
74916413a6aSBarry Smith   ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)ls),((PetscObject) ls)->options,((PetscObject)ls)->prefix,name,&viewer,&format,&flg);CHKERRQ(ierr);
750dcf2fd19SBarry Smith   if (flg) {
751d12e167eSBarry Smith     PetscViewerAndFormat *vf;
752d12e167eSBarry Smith     ierr = PetscViewerAndFormatCreate(viewer,format,&vf);CHKERRQ(ierr);
753d12e167eSBarry Smith     ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr);
754dcf2fd19SBarry Smith     if (monitorsetup) {
755d12e167eSBarry Smith       ierr = (*monitorsetup)(ls,vf);CHKERRQ(ierr);
756dcf2fd19SBarry Smith     }
757d12e167eSBarry Smith     ierr = SNESLineSearchMonitorSet(ls,(PetscErrorCode (*)(SNESLineSearch,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);CHKERRQ(ierr);
758dcf2fd19SBarry Smith   }
759dcf2fd19SBarry Smith   PetscFunctionReturn(0);
760dcf2fd19SBarry Smith }
761dcf2fd19SBarry Smith 
762f40b411bSPeter Brune /*@
763f1c6b773SPeter Brune    SNESLineSearchSetFromOptions - Sets options for the line search
764f40b411bSPeter Brune 
765f40b411bSPeter Brune    Input Parameters:
7668e557f58SPeter Brune .  linesearch - linesearch context
767f40b411bSPeter Brune 
768f40b411bSPeter Brune    Options Database Keys:
769d4c6564cSPatrick Farrell + -snes_linesearch_type <type> - basic, bt, l2, cp, nleqerr, shell
7703c7d6663SPeter Brune . -snes_linesearch_order <order> - 1, 2, 3.  Most types only support certain orders (bt supports 2 or 3)
7711fe24845SBarry Smith . -snes_linesearch_norms   - Turn on/off the linesearch norms for the basic linesearch typem (SNESLineSearchSetComputeNorms())
77271eef1aeSPeter Brune . -snes_linesearch_minlambda - The minimum step length
7731a9b3a06SPeter Brune . -snes_linesearch_maxstep - The maximum step size
7741a9b3a06SPeter Brune . -snes_linesearch_rtol - Relative tolerance for iterative line searches
7751a9b3a06SPeter Brune . -snes_linesearch_atol - Absolute tolerance for iterative line searches
7761a9b3a06SPeter Brune . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches
7771a9b3a06SPeter Brune . -snes_linesearch_max_it - The number of iterations for iterative line searches
778dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches
779dcf2fd19SBarry Smith . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine
7808e557f58SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter
781cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess.
7821a9b3a06SPeter Brune . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method
783d8d34be6SBarry Smith - -snes_linesearch_precheck_picard_angle - Angle used in Picard precheck method
784f40b411bSPeter Brune 
785f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
786f40b411bSPeter Brune 
787f40b411bSPeter Brune    Level: intermediate
788f40b411bSPeter Brune 
7891fe24845SBarry Smith .seealso: SNESLineSearchCreate(), SNESLineSearchSetOrder(), SNESLineSearchSetType(), SNESLineSearchSetTolerances(), SNESLineSearchSetDamping(), SNESLineSearchPreCheckPicard(),
7901fe24845SBarry Smith           SNESLineSearchType, SNESLineSearchSetComputeNorms()
791f40b411bSPeter Brune @*/
792bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch)
793bf388a1fSBarry Smith {
794bf7f4e0aSPeter Brune   PetscErrorCode    ierr;
7951a4f838cSPeter Brune   const char        *deft = SNESLINESEARCHBASIC;
796bf7f4e0aSPeter Brune   char              type[256];
797bf7f4e0aSPeter Brune   PetscBool         flg, set;
798dcf2fd19SBarry Smith   PetscViewer       viewer;
799bf388a1fSBarry Smith 
800bf7f4e0aSPeter Brune   PetscFunctionBegin;
8010f51fdf8SToby Isaac   ierr = SNESLineSearchRegisterAll();CHKERRQ(ierr);
802bf7f4e0aSPeter Brune 
803bf7f4e0aSPeter Brune   ierr = PetscObjectOptionsBegin((PetscObject)linesearch);CHKERRQ(ierr);
804f5af7f23SKarl Rupp   if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name;
805a264d7a6SBarry Smith   ierr = PetscOptionsFList("-snes_linesearch_type","Linesearch type","SNESLineSearchSetType",SNESLineSearchList,deft,type,256,&flg);CHKERRQ(ierr);
806bf7f4e0aSPeter Brune   if (flg) {
807f1c6b773SPeter Brune     ierr = SNESLineSearchSetType(linesearch,type);CHKERRQ(ierr);
808bf7f4e0aSPeter Brune   } else if (!((PetscObject)linesearch)->type_name) {
809f1c6b773SPeter Brune     ierr = SNESLineSearchSetType(linesearch,deft);CHKERRQ(ierr);
810bf7f4e0aSPeter Brune   }
811bf7f4e0aSPeter Brune 
81216413a6aSBarry Smith   ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)linesearch),((PetscObject) linesearch)->options,((PetscObject)linesearch)->prefix,"-snes_linesearch_monitor",&viewer,NULL,&set);CHKERRQ(ierr);
813dcf2fd19SBarry Smith   if (set) {
814dcf2fd19SBarry Smith     ierr = SNESLineSearchSetDefaultMonitor(linesearch,viewer);CHKERRQ(ierr);
815dcf2fd19SBarry Smith     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
816dcf2fd19SBarry Smith   }
817dcf2fd19SBarry Smith   ierr = SNESLineSearchMonitorSetFromOptions(linesearch,"-snes_linesearch_monitor_solution_update","View correction at each iteration","SNESLineSearchMonitorSolutionUpdate",SNESLineSearchMonitorSolutionUpdate,NULL);CHKERRQ(ierr);
818bf7f4e0aSPeter Brune 
8191a9b3a06SPeter Brune   /* tolerances */
82094ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_minlambda","Minimum step length","SNESLineSearchSetTolerances",linesearch->steptol,&linesearch->steptol,NULL);CHKERRQ(ierr);
82194ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_maxstep","Maximum step size","SNESLineSearchSetTolerances",linesearch->maxstep,&linesearch->maxstep,NULL);CHKERRQ(ierr);
82294ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_rtol","Relative tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->rtol,&linesearch->rtol,NULL);CHKERRQ(ierr);
82394ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_atol","Absolute tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->atol,&linesearch->atol,NULL);CHKERRQ(ierr);
82494ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_ltol","Change in lambda tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->ltol,&linesearch->ltol,NULL);CHKERRQ(ierr);
82594ae4db5SBarry Smith   ierr = PetscOptionsInt("-snes_linesearch_max_it","Maximum iterations for iterative line searches","SNESLineSearchSetTolerances",linesearch->max_its,&linesearch->max_its,NULL);CHKERRQ(ierr);
8267a35526eSPeter Brune 
8271a9b3a06SPeter Brune   /* damping parameters */
82894ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_damping","Line search damping and initial step guess","SNESLineSearchSetDamping",linesearch->damping,&linesearch->damping,NULL);CHKERRQ(ierr);
8291a9b3a06SPeter Brune 
83094ae4db5SBarry Smith   ierr = PetscOptionsBool("-snes_linesearch_keeplambda","Use previous lambda as damping","SNESLineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,NULL);CHKERRQ(ierr);
8311a9b3a06SPeter Brune 
8321a9b3a06SPeter Brune   /* precheck */
8337a35526eSPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_precheck_picard","Use a correction that sometimes improves convergence of Picard iteration","SNESLineSearchPreCheckPicard",flg,&flg,&set);CHKERRQ(ierr);
8347a35526eSPeter Brune   if (set) {
8357a35526eSPeter Brune     if (flg) {
8367a35526eSPeter Brune       linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */
837f5af7f23SKarl Rupp 
8387a35526eSPeter Brune       ierr = PetscOptionsReal("-snes_linesearch_precheck_picard_angle","Maximum angle at which to activate the correction",
8390298fd71SBarry Smith                               "none",linesearch->precheck_picard_angle,&linesearch->precheck_picard_angle,NULL);CHKERRQ(ierr);
8407a35526eSPeter Brune       ierr = SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&linesearch->precheck_picard_angle);CHKERRQ(ierr);
8417a35526eSPeter Brune     } else {
8420298fd71SBarry Smith       ierr = SNESLineSearchSetPreCheck(linesearch,NULL,NULL);CHKERRQ(ierr);
8437a35526eSPeter Brune     }
8447a35526eSPeter Brune   }
84594ae4db5SBarry Smith   ierr = PetscOptionsInt("-snes_linesearch_order","Order of approximation used in the line search","SNESLineSearchSetOrder",linesearch->order,&linesearch->order,NULL);CHKERRQ(ierr);
84694ae4db5SBarry Smith   ierr = PetscOptionsBool("-snes_linesearch_norms","Compute final norms in line search","SNESLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,NULL);CHKERRQ(ierr);
8477a35526eSPeter Brune 
8485a9b6599SPeter Brune   if (linesearch->ops->setfromoptions) {
849e55864a3SBarry Smith     (*linesearch->ops->setfromoptions)(PetscOptionsObject,linesearch);CHKERRQ(ierr);
8505a9b6599SPeter Brune   }
8515a9b6599SPeter Brune 
8520633abcbSJed Brown   ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)linesearch);CHKERRQ(ierr);
853bf7f4e0aSPeter Brune   ierr = PetscOptionsEnd();CHKERRQ(ierr);
854bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
855bf7f4e0aSPeter Brune }
856bf7f4e0aSPeter Brune 
857f40b411bSPeter Brune /*@
858f190f2fcSBarry Smith    SNESLineSearchView - Prints useful information about the line search
859f40b411bSPeter Brune 
860f40b411bSPeter Brune    Input Parameters:
8618e557f58SPeter Brune .  linesearch - linesearch context
862f40b411bSPeter Brune 
863f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
864f40b411bSPeter Brune 
865f40b411bSPeter Brune    Level: intermediate
866f40b411bSPeter Brune 
867dcf2fd19SBarry Smith .seealso: SNESLineSearchCreate()
868f40b411bSPeter Brune @*/
869bf388a1fSBarry Smith PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer)
870bf388a1fSBarry Smith {
8717f1410a3SPeter Brune   PetscErrorCode ierr;
8727f1410a3SPeter Brune   PetscBool      iascii;
873bf388a1fSBarry Smith 
874bf7f4e0aSPeter Brune   PetscFunctionBegin;
8757f1410a3SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
8767f1410a3SPeter Brune   if (!viewer) {
877ce94432eSBarry Smith     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch),&viewer);CHKERRQ(ierr);
8787f1410a3SPeter Brune   }
8797f1410a3SPeter Brune   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
8807f1410a3SPeter Brune   PetscCheckSameComm(linesearch,1,viewer,2);
881f40b411bSPeter Brune 
882251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
8837f1410a3SPeter Brune   if (iascii) {
884dae58748SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)linesearch,viewer);CHKERRQ(ierr);
8857f1410a3SPeter Brune     if (linesearch->ops->view) {
8867f1410a3SPeter Brune       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
8877f1410a3SPeter Brune       ierr = (*linesearch->ops->view)(linesearch,viewer);CHKERRQ(ierr);
8887f1410a3SPeter Brune       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
8897f1410a3SPeter Brune     }
890c69d1a72SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep,(double)linesearch->steptol);CHKERRQ(ierr);
891c69d1a72SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol,(double)linesearch->atol,(double)linesearch->ltol);CHKERRQ(ierr);
8927f1410a3SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D\n", linesearch->max_its);CHKERRQ(ierr);
8936b2b7091SBarry Smith     if (linesearch->ops->precheck) {
8946b2b7091SBarry Smith       if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) {
8957f1410a3SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  using precheck step to speed up Picard convergence\n", linesearch->max_its);CHKERRQ(ierr);
8967f1410a3SPeter Brune       } else {
8977f1410a3SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  using user-defined precheck step\n", linesearch->max_its);CHKERRQ(ierr);
8987f1410a3SPeter Brune       }
8997f1410a3SPeter Brune     }
9006b2b7091SBarry Smith     if (linesearch->ops->postcheck) {
9017f1410a3SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer,"  using user-defined postcheck step\n", linesearch->max_its);CHKERRQ(ierr);
9027f1410a3SPeter Brune     }
9037f1410a3SPeter Brune   }
904bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
905bf7f4e0aSPeter Brune }
906bf7f4e0aSPeter Brune 
907ea5d4fccSPeter Brune /*@C
908f1c6b773SPeter Brune    SNESLineSearchSetType - Sets the linesearch type
909f40b411bSPeter Brune 
910f190f2fcSBarry Smith    Logically Collective on SNESLineSearch
911f190f2fcSBarry Smith 
912f40b411bSPeter Brune    Input Parameters:
9138e557f58SPeter Brune +  linesearch - linesearch context
914f40b411bSPeter Brune -  type - The type of line search to be used
915f40b411bSPeter Brune 
916cd7522eaSPeter Brune    Available Types:
9171fe24845SBarry Smith +  SNESLINESEARCHBASIC - Simple damping line search, defaults to using the full Newton step
9181fe24845SBarry Smith .  SNESLINESEARCHBT - Backtracking line search over the L2 norm of the function
9191fe24845SBarry Smith .  SNESLINESEARCHL2 - Secant line search over the L2 norm of the function
9201fe24845SBarry Smith .  SNESLINESEARCHCP - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x)
9211fe24845SBarry Smith .  SNESLINESEARCHNLEQERR - Affine-covariant error-oriented linesearch
9221fe24845SBarry Smith -  SNESLINESEARCHSHELL - User provided SNESLineSearch implementation
9231fe24845SBarry Smith 
9241fe24845SBarry Smith    Options Database:
9251fe24845SBarry Smith .  -snes_linesearch_type <type> - basic, bt, l2, cp, nleqerr, shell
926cd7522eaSPeter Brune 
927f40b411bSPeter Brune    Level: intermediate
928f40b411bSPeter Brune 
9291fe24845SBarry Smith .seealso: SNESLineSearchCreate(), SNESLineSearchType, SNESLineSearchSetFromOptions()
930f40b411bSPeter Brune @*/
93119fd82e9SBarry Smith PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type)
932bf7f4e0aSPeter Brune {
933f1c6b773SPeter Brune   PetscErrorCode ierr,(*r)(SNESLineSearch);
934bf7f4e0aSPeter Brune   PetscBool      match;
935bf7f4e0aSPeter Brune 
936bf7f4e0aSPeter Brune   PetscFunctionBegin;
937f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
938bf7f4e0aSPeter Brune   PetscValidCharPointer(type,2);
939bf7f4e0aSPeter Brune 
940251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)linesearch,type,&match);CHKERRQ(ierr);
941bf7f4e0aSPeter Brune   if (match) PetscFunctionReturn(0);
942bf7f4e0aSPeter Brune 
9431c9cd337SJed Brown   ierr = PetscFunctionListFind(SNESLineSearchList,type,&r);CHKERRQ(ierr);
944bf7f4e0aSPeter Brune   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type);
945bf7f4e0aSPeter Brune   /* Destroy the previous private linesearch context */
946bf7f4e0aSPeter Brune   if (linesearch->ops->destroy) {
947bf7f4e0aSPeter Brune     ierr = (*(linesearch)->ops->destroy)(linesearch);CHKERRQ(ierr);
948f5af7f23SKarl Rupp 
9490298fd71SBarry Smith     linesearch->ops->destroy = NULL;
950bf7f4e0aSPeter Brune   }
951f1c6b773SPeter Brune   /* Reinitialize function pointers in SNESLineSearchOps structure */
952bf7f4e0aSPeter Brune   linesearch->ops->apply          = 0;
953bf7f4e0aSPeter Brune   linesearch->ops->view           = 0;
954bf7f4e0aSPeter Brune   linesearch->ops->setfromoptions = 0;
955bf7f4e0aSPeter Brune   linesearch->ops->destroy        = 0;
956bf7f4e0aSPeter Brune 
957bf7f4e0aSPeter Brune   ierr = PetscObjectChangeTypeName((PetscObject)linesearch,type);CHKERRQ(ierr);
958bf7f4e0aSPeter Brune   ierr = (*r)(linesearch);CHKERRQ(ierr);
959bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
960bf7f4e0aSPeter Brune }
961bf7f4e0aSPeter Brune 
962f40b411bSPeter Brune /*@
96378bcb3b5SPeter Brune    SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation.
964f40b411bSPeter Brune 
965f40b411bSPeter Brune    Input Parameters:
9668e557f58SPeter Brune +  linesearch - linesearch context
967f40b411bSPeter Brune -  snes - The snes instance
968f40b411bSPeter Brune 
96978bcb3b5SPeter Brune    Level: developer
97078bcb3b5SPeter Brune 
97178bcb3b5SPeter Brune    Notes:
972f190f2fcSBarry Smith    This happens automatically when the line search is obtained/created with
9737601faf0SJed Brown    SNESGetLineSearch().  This routine is therefore mainly called within SNES
97478bcb3b5SPeter Brune    implementations.
975f40b411bSPeter Brune 
9768141a3b9SPeter Brune    Level: developer
977f40b411bSPeter Brune 
978cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
979f40b411bSPeter Brune @*/
980bf388a1fSBarry Smith PetscErrorCode  SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes)
981bf388a1fSBarry Smith {
982bf7f4e0aSPeter Brune   PetscFunctionBegin;
983f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
984bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,2);
985bf7f4e0aSPeter Brune   linesearch->snes = snes;
986bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
987bf7f4e0aSPeter Brune }
988bf7f4e0aSPeter Brune 
989f40b411bSPeter Brune /*@
9908141a3b9SPeter Brune    SNESLineSearchGetSNES - Gets the SNES instance associated with the line search.
9918141a3b9SPeter Brune    Having an associated SNES is necessary because most line search implementations must be able to
9928141a3b9SPeter Brune    evaluate the function using SNESComputeFunction() for the associated SNES.  This routine
9938141a3b9SPeter Brune    is used in the line search implementations when one must get this associated SNES instance.
994f40b411bSPeter Brune 
995f40b411bSPeter Brune    Input Parameters:
9968e557f58SPeter Brune .  linesearch - linesearch context
997f40b411bSPeter Brune 
998f40b411bSPeter Brune    Output Parameters:
999f40b411bSPeter Brune .  snes - The snes instance
1000f40b411bSPeter Brune 
10018141a3b9SPeter Brune    Level: developer
1002f40b411bSPeter Brune 
1003cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
1004f40b411bSPeter Brune @*/
1005bf388a1fSBarry Smith PetscErrorCode  SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes)
1006bf388a1fSBarry Smith {
1007bf7f4e0aSPeter Brune   PetscFunctionBegin;
1008f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
10096a388c36SPeter Brune   PetscValidPointer(snes, 2);
1010bf7f4e0aSPeter Brune   *snes = linesearch->snes;
1011bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1012bf7f4e0aSPeter Brune }
1013bf7f4e0aSPeter Brune 
1014f40b411bSPeter Brune /*@
1015f1c6b773SPeter Brune    SNESLineSearchGetLambda - Gets the last linesearch steplength discovered.
1016f40b411bSPeter Brune 
1017f40b411bSPeter Brune    Input Parameters:
10188e557f58SPeter Brune .  linesearch - linesearch context
1019f40b411bSPeter Brune 
1020f40b411bSPeter Brune    Output Parameters:
1021cd7522eaSPeter Brune .  lambda - The last steplength computed during SNESLineSearchApply()
1022f40b411bSPeter Brune 
102378bcb3b5SPeter Brune    Level: advanced
102478bcb3b5SPeter Brune 
10258e557f58SPeter Brune    Notes:
10268e557f58SPeter Brune    This is useful in methods where the solver is ill-scaled and
102778bcb3b5SPeter Brune    requires some adaptive notion of the difference in scale between the
102878bcb3b5SPeter Brune    solution and the function.  For instance, SNESQN may be scaled by the
102978bcb3b5SPeter Brune    line search lambda using the argument -snes_qn_scaling ls.
103078bcb3b5SPeter Brune 
1031cd7522eaSPeter Brune .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetDamping(), SNESLineSearchApply()
1032f40b411bSPeter Brune @*/
1033f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda)
10346a388c36SPeter Brune {
10356a388c36SPeter Brune   PetscFunctionBegin;
1036f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1037*534a8f05SLisandro Dalcin   PetscValidRealPointer(lambda, 2);
10386a388c36SPeter Brune   *lambda = linesearch->lambda;
10396a388c36SPeter Brune   PetscFunctionReturn(0);
10406a388c36SPeter Brune }
10416a388c36SPeter Brune 
1042f40b411bSPeter Brune /*@
1043f1c6b773SPeter Brune    SNESLineSearchSetLambda - Sets the linesearch steplength.
1044f40b411bSPeter Brune 
1045f40b411bSPeter Brune    Input Parameters:
10468e557f58SPeter Brune +  linesearch - linesearch context
1047f40b411bSPeter Brune -  lambda - The last steplength.
1048f40b411bSPeter Brune 
1049cd7522eaSPeter Brune    Notes:
1050f190f2fcSBarry Smith    This routine is typically used within implementations of SNESLineSearchApply()
1051cd7522eaSPeter Brune    to set the final steplength.  This routine (and SNESLineSearchGetLambda()) were
1052cd7522eaSPeter Brune    added in order to facilitate Quasi-Newton methods that use the previous steplength
1053cd7522eaSPeter Brune    as an inner scaling parameter.
1054cd7522eaSPeter Brune 
105578bcb3b5SPeter Brune    Level: advanced
1056f40b411bSPeter Brune 
1057f1c6b773SPeter Brune .seealso: SNESLineSearchGetLambda()
1058f40b411bSPeter Brune @*/
1059f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
10606a388c36SPeter Brune {
10616a388c36SPeter Brune   PetscFunctionBegin;
1062f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
10636a388c36SPeter Brune   linesearch->lambda = lambda;
10646a388c36SPeter Brune   PetscFunctionReturn(0);
10656a388c36SPeter Brune }
10666a388c36SPeter Brune 
1067f40b411bSPeter Brune /*@
10683c7d6663SPeter Brune    SNESLineSearchGetTolerances - Gets the tolerances for the linesearch.  These include
106978bcb3b5SPeter Brune    tolerances for the relative and absolute change in the function norm, the change
107078bcb3b5SPeter Brune    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
107178bcb3b5SPeter Brune    and the maximum number of iterations the line search procedure may take.
1072f40b411bSPeter Brune 
1073f40b411bSPeter Brune    Input Parameters:
10748e557f58SPeter Brune .  linesearch - linesearch context
1075f40b411bSPeter Brune 
1076f40b411bSPeter Brune    Output Parameters:
1077516fe3c3SPeter Brune +  steptol - The minimum steplength
10786cc8e53bSPeter Brune .  maxstep - The maximum steplength
1079516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
1080516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
1081516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
1082516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
1083f40b411bSPeter Brune 
108478bcb3b5SPeter Brune    Level: intermediate
108578bcb3b5SPeter Brune 
108678bcb3b5SPeter Brune    Notes:
108778bcb3b5SPeter Brune    Different line searches may implement these parameters slightly differently as
10883c7d6663SPeter Brune    the type requires.
1089516fe3c3SPeter Brune 
1090f1c6b773SPeter Brune .seealso: SNESLineSearchSetTolerances()
1091f40b411bSPeter Brune @*/
1092f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
10936a388c36SPeter Brune {
10946a388c36SPeter Brune   PetscFunctionBegin;
1095f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1096516fe3c3SPeter Brune   if (steptol) {
1097*534a8f05SLisandro Dalcin     PetscValidRealPointer(steptol, 2);
10986a388c36SPeter Brune     *steptol = linesearch->steptol;
1099516fe3c3SPeter Brune   }
1100516fe3c3SPeter Brune   if (maxstep) {
1101*534a8f05SLisandro Dalcin     PetscValidRealPointer(maxstep, 3);
1102516fe3c3SPeter Brune     *maxstep = linesearch->maxstep;
1103516fe3c3SPeter Brune   }
1104516fe3c3SPeter Brune   if (rtol) {
1105*534a8f05SLisandro Dalcin     PetscValidRealPointer(rtol, 4);
1106516fe3c3SPeter Brune     *rtol = linesearch->rtol;
1107516fe3c3SPeter Brune   }
1108516fe3c3SPeter Brune   if (atol) {
1109*534a8f05SLisandro Dalcin     PetscValidRealPointer(atol, 5);
1110516fe3c3SPeter Brune     *atol = linesearch->atol;
1111516fe3c3SPeter Brune   }
1112516fe3c3SPeter Brune   if (ltol) {
1113*534a8f05SLisandro Dalcin     PetscValidRealPointer(ltol, 6);
1114516fe3c3SPeter Brune     *ltol = linesearch->ltol;
1115516fe3c3SPeter Brune   }
1116516fe3c3SPeter Brune   if (max_its) {
1117*534a8f05SLisandro Dalcin     PetscValidIntPointer(max_its, 7);
1118516fe3c3SPeter Brune     *max_its = linesearch->max_its;
1119516fe3c3SPeter Brune   }
11206a388c36SPeter Brune   PetscFunctionReturn(0);
11216a388c36SPeter Brune }
11226a388c36SPeter Brune 
1123f40b411bSPeter Brune /*@
11243c7d6663SPeter Brune    SNESLineSearchSetTolerances -  Gets the tolerances for the linesearch.  These include
112578bcb3b5SPeter Brune    tolerances for the relative and absolute change in the function norm, the change
112678bcb3b5SPeter Brune    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
112778bcb3b5SPeter Brune    and the maximum number of iterations the line search procedure may take.
1128f40b411bSPeter Brune 
1129f40b411bSPeter Brune    Input Parameters:
11308e557f58SPeter Brune +  linesearch - linesearch context
1131516fe3c3SPeter Brune .  steptol - The minimum steplength
11326cc8e53bSPeter Brune .  maxstep - The maximum steplength
1133516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
1134516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
1135516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
1136516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
1137f40b411bSPeter Brune 
113878bcb3b5SPeter Brune    Notes:
11393c7d6663SPeter Brune    The user may choose to not set any of the tolerances using PETSC_DEFAULT in
114078bcb3b5SPeter Brune    place of an argument.
1141f40b411bSPeter Brune 
114278bcb3b5SPeter Brune    Level: intermediate
1143516fe3c3SPeter Brune 
1144f1c6b773SPeter Brune .seealso: SNESLineSearchGetTolerances()
1145f40b411bSPeter Brune @*/
1146f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its)
11476a388c36SPeter Brune {
11486a388c36SPeter Brune   PetscFunctionBegin;
1149f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1150d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,steptol,2);
1151d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,maxstep,3);
1152d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,rtol,4);
1153d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,atol,5);
1154d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,ltol,6);
1155d3952184SSatish Balay   PetscValidLogicalCollectiveInt(linesearch,max_its,7);
1156d3952184SSatish Balay 
1157d3952184SSatish Balay   if (steptol!= PETSC_DEFAULT) {
1158ce94432eSBarry Smith     if (steptol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Minimum step length %14.12e must be non-negative",(double)steptol);
11596a388c36SPeter Brune     linesearch->steptol = steptol;
1160d3952184SSatish Balay   }
1161d3952184SSatish Balay 
1162d3952184SSatish Balay   if (maxstep!= PETSC_DEFAULT) {
1163ce94432eSBarry Smith     if (maxstep < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum step length %14.12e must be non-negative",(double)maxstep);
1164516fe3c3SPeter Brune     linesearch->maxstep = maxstep;
1165d3952184SSatish Balay   }
1166d3952184SSatish Balay 
1167d3952184SSatish Balay   if (rtol != PETSC_DEFAULT) {
1168ce94432eSBarry 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);
1169516fe3c3SPeter Brune     linesearch->rtol = rtol;
1170d3952184SSatish Balay   }
1171d3952184SSatish Balay 
1172d3952184SSatish Balay   if (atol != PETSC_DEFAULT) {
1173ce94432eSBarry Smith     if (atol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %14.12e must be non-negative",(double)atol);
1174516fe3c3SPeter Brune     linesearch->atol = atol;
1175d3952184SSatish Balay   }
1176d3952184SSatish Balay 
1177d3952184SSatish Balay   if (ltol != PETSC_DEFAULT) {
1178ce94432eSBarry Smith     if (ltol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Labmda tolerance %14.12e must be non-negative",(double)ltol);
1179516fe3c3SPeter Brune     linesearch->ltol = ltol;
1180d3952184SSatish Balay   }
1181d3952184SSatish Balay 
1182d3952184SSatish Balay   if (max_its != PETSC_DEFAULT) {
1183ce94432eSBarry Smith     if (max_its < 0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",max_its);
1184516fe3c3SPeter Brune     linesearch->max_its = max_its;
1185d3952184SSatish Balay   }
11866a388c36SPeter Brune   PetscFunctionReturn(0);
11876a388c36SPeter Brune }
11886a388c36SPeter Brune 
1189f40b411bSPeter Brune /*@
1190f1c6b773SPeter Brune    SNESLineSearchGetDamping - Gets the line search damping parameter.
1191f40b411bSPeter Brune 
1192f40b411bSPeter Brune    Input Parameters:
11938e557f58SPeter Brune .  linesearch - linesearch context
1194f40b411bSPeter Brune 
1195f40b411bSPeter Brune    Output Parameters:
11968e557f58SPeter Brune .  damping - The damping parameter
1197f40b411bSPeter Brune 
119878bcb3b5SPeter Brune    Level: advanced
1199f40b411bSPeter Brune 
120078bcb3b5SPeter Brune .seealso: SNESLineSearchGetStepTolerance(), SNESQN
1201f40b411bSPeter Brune @*/
1202f40b411bSPeter Brune 
1203f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping)
12046a388c36SPeter Brune {
12056a388c36SPeter Brune   PetscFunctionBegin;
1206f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1207*534a8f05SLisandro Dalcin   PetscValidRealPointer(damping, 2);
12086a388c36SPeter Brune   *damping = linesearch->damping;
12096a388c36SPeter Brune   PetscFunctionReturn(0);
12106a388c36SPeter Brune }
12116a388c36SPeter Brune 
1212f40b411bSPeter Brune /*@
1213f1c6b773SPeter Brune    SNESLineSearchSetDamping - Sets the line search damping paramter.
1214f40b411bSPeter Brune 
1215f40b411bSPeter Brune    Input Parameters:
121603fd4120SBarry Smith +  linesearch - linesearch context
121703fd4120SBarry Smith -  damping - The damping parameter
1218f40b411bSPeter Brune 
121903fd4120SBarry Smith    Options Database:
122003fd4120SBarry Smith .   -snes_linesearch_damping
1221f40b411bSPeter Brune    Level: intermediate
1222f40b411bSPeter Brune 
1223cd7522eaSPeter Brune    Notes:
1224cd7522eaSPeter Brune    The basic line search merely takes the update step scaled by the damping parameter.
1225cd7522eaSPeter Brune    The use of the damping parameter in the l2 and cp line searches is much more subtle;
122678bcb3b5SPeter Brune    it is used as a starting point in calculating the secant step. However, the eventual
1227cd7522eaSPeter Brune    step may be of greater length than the damping parameter.  In the bt line search it is
1228cd7522eaSPeter Brune    used as the maximum possible step length, as the bt line search only backtracks.
1229cd7522eaSPeter Brune 
1230f1c6b773SPeter Brune .seealso: SNESLineSearchGetDamping()
1231f40b411bSPeter Brune @*/
1232f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping)
12336a388c36SPeter Brune {
12346a388c36SPeter Brune   PetscFunctionBegin;
1235f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
12366a388c36SPeter Brune   linesearch->damping = damping;
12376a388c36SPeter Brune   PetscFunctionReturn(0);
12386a388c36SPeter Brune }
12396a388c36SPeter Brune 
124059405d5eSPeter Brune /*@
124159405d5eSPeter Brune    SNESLineSearchGetOrder - Gets the line search approximation order.
124259405d5eSPeter Brune 
124359405d5eSPeter Brune    Input Parameters:
124478bcb3b5SPeter Brune .  linesearch - linesearch context
124559405d5eSPeter Brune 
124659405d5eSPeter Brune    Output Parameters:
12478e557f58SPeter Brune .  order - The order
124859405d5eSPeter Brune 
124978bcb3b5SPeter Brune    Possible Values for order:
12503c7d6663SPeter Brune +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
12513c7d6663SPeter Brune .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
12523c7d6663SPeter Brune -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
125378bcb3b5SPeter Brune 
125459405d5eSPeter Brune    Level: intermediate
125559405d5eSPeter Brune 
125659405d5eSPeter Brune .seealso: SNESLineSearchSetOrder()
125759405d5eSPeter Brune @*/
125859405d5eSPeter Brune 
1259b000cd8dSPeter Brune PetscErrorCode  SNESLineSearchGetOrder(SNESLineSearch linesearch,PetscInt *order)
126059405d5eSPeter Brune {
126159405d5eSPeter Brune   PetscFunctionBegin;
126259405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1263*534a8f05SLisandro Dalcin   PetscValidIntPointer(order, 2);
126459405d5eSPeter Brune   *order = linesearch->order;
126559405d5eSPeter Brune   PetscFunctionReturn(0);
126659405d5eSPeter Brune }
126759405d5eSPeter Brune 
126859405d5eSPeter Brune /*@
12691f8196a2SJed Brown    SNESLineSearchSetOrder - Sets the maximum order of the polynomial fit used in the line search
127059405d5eSPeter Brune 
127159405d5eSPeter Brune    Input Parameters:
1272a2b725a8SWilliam Gropp +  linesearch - linesearch context
1273a2b725a8SWilliam Gropp -  order - The damping parameter
127459405d5eSPeter Brune 
127559405d5eSPeter Brune    Level: intermediate
127659405d5eSPeter Brune 
127778bcb3b5SPeter Brune    Possible Values for order:
12783c7d6663SPeter Brune +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
12793c7d6663SPeter Brune .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
12803c7d6663SPeter Brune -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
128178bcb3b5SPeter Brune 
128259405d5eSPeter Brune    Notes:
128359405d5eSPeter Brune    Variable orders are supported by the following line searches:
128478bcb3b5SPeter Brune +  bt - cubic and quadratic
128578bcb3b5SPeter Brune -  cp - linear and quadratic
128659405d5eSPeter Brune 
12871f8196a2SJed Brown .seealso: SNESLineSearchGetOrder(), SNESLineSearchSetDamping()
128859405d5eSPeter Brune @*/
1289b000cd8dSPeter Brune PetscErrorCode  SNESLineSearchSetOrder(SNESLineSearch linesearch,PetscInt order)
129059405d5eSPeter Brune {
129159405d5eSPeter Brune   PetscFunctionBegin;
129259405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
129359405d5eSPeter Brune   linesearch->order = order;
129459405d5eSPeter Brune   PetscFunctionReturn(0);
129559405d5eSPeter Brune }
129659405d5eSPeter Brune 
1297f40b411bSPeter Brune /*@
1298f1c6b773SPeter Brune    SNESLineSearchGetNorms - Gets the norms for for X, Y, and F.
1299f40b411bSPeter Brune 
1300f40b411bSPeter Brune    Input Parameters:
130178bcb3b5SPeter Brune .  linesearch - linesearch context
1302f40b411bSPeter Brune 
1303f40b411bSPeter Brune    Output Parameters:
1304f40b411bSPeter Brune +  xnorm - The norm of the current solution
1305f40b411bSPeter Brune .  fnorm - The norm of the current function
1306f40b411bSPeter Brune -  ynorm - The norm of the current update
1307f40b411bSPeter Brune 
1308cd7522eaSPeter Brune    Notes:
1309cd7522eaSPeter Brune    This function is mainly called from SNES implementations.
1310cd7522eaSPeter Brune 
131178bcb3b5SPeter Brune    Level: developer
1312f40b411bSPeter Brune 
1313f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms() SNESLineSearchGetVecs()
1314f40b411bSPeter Brune @*/
1315f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm)
1316bf7f4e0aSPeter Brune {
1317bf7f4e0aSPeter Brune   PetscFunctionBegin;
1318f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1319f5af7f23SKarl Rupp   if (xnorm) *xnorm = linesearch->xnorm;
1320f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
1321f5af7f23SKarl Rupp   if (ynorm) *ynorm = linesearch->ynorm;
1322bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1323bf7f4e0aSPeter Brune }
1324bf7f4e0aSPeter Brune 
1325f40b411bSPeter Brune /*@
1326f1c6b773SPeter Brune    SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F.
1327f40b411bSPeter Brune 
1328f40b411bSPeter Brune    Input Parameters:
132978bcb3b5SPeter Brune +  linesearch - linesearch context
1330f40b411bSPeter Brune .  xnorm - The norm of the current solution
1331f40b411bSPeter Brune .  fnorm - The norm of the current function
1332f40b411bSPeter Brune -  ynorm - The norm of the current update
1333f40b411bSPeter Brune 
133478bcb3b5SPeter Brune    Level: advanced
1335f40b411bSPeter Brune 
1336f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1337f40b411bSPeter Brune @*/
1338f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
13396a388c36SPeter Brune {
13406a388c36SPeter Brune   PetscFunctionBegin;
1341f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
13426a388c36SPeter Brune   linesearch->xnorm = xnorm;
13436a388c36SPeter Brune   linesearch->fnorm = fnorm;
13446a388c36SPeter Brune   linesearch->ynorm = ynorm;
13456a388c36SPeter Brune   PetscFunctionReturn(0);
13466a388c36SPeter Brune }
13476a388c36SPeter Brune 
1348f40b411bSPeter Brune /*@
1349f1c6b773SPeter Brune    SNESLineSearchComputeNorms - Computes the norms of X, F, and Y.
1350f40b411bSPeter Brune 
1351f40b411bSPeter Brune    Input Parameters:
135278bcb3b5SPeter Brune .  linesearch - linesearch context
1353f40b411bSPeter Brune 
1354f40b411bSPeter Brune    Options Database Keys:
13558e557f58SPeter Brune .   -snes_linesearch_norms - turn norm computation on or off
1356f40b411bSPeter Brune 
1357f40b411bSPeter Brune    Level: intermediate
1358f40b411bSPeter Brune 
135978bcb3b5SPeter Brune .seealso: SNESLineSearchGetNorms, SNESLineSearchSetNorms(), SNESLineSearchSetComputeNorms()
1360f40b411bSPeter Brune @*/
1361f1c6b773SPeter Brune PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
13626a388c36SPeter Brune {
13636a388c36SPeter Brune   PetscErrorCode ierr;
13649bd66eb0SPeter Brune   SNES           snes;
1365bf388a1fSBarry Smith 
13666a388c36SPeter Brune   PetscFunctionBegin;
13676a388c36SPeter Brune   if (linesearch->norms) {
13689bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
1369f1c6b773SPeter Brune       ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
13709bd66eb0SPeter Brune       ierr = VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
13719bd66eb0SPeter Brune       ierr = VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
13729bd66eb0SPeter Brune       ierr = (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);CHKERRQ(ierr);
13739bd66eb0SPeter Brune     } else {
13746a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_func,   NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
13756a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_sol,    NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
13766a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
13776a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_func,     NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
13786a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_sol,      NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
13796a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_update,   NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
13806a388c36SPeter Brune     }
13819bd66eb0SPeter Brune   }
13826a388c36SPeter Brune   PetscFunctionReturn(0);
13836a388c36SPeter Brune }
13846a388c36SPeter Brune 
13856f263ca3SPeter Brune /*@
13866f263ca3SPeter Brune    SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search.
13876f263ca3SPeter Brune 
13886f263ca3SPeter Brune    Input Parameters:
138978bcb3b5SPeter Brune +  linesearch  - linesearch context
139078bcb3b5SPeter Brune -  flg  - indicates whether or not to compute norms
13916f263ca3SPeter Brune 
13926f263ca3SPeter Brune    Options Database Keys:
13931f8196a2SJed Brown .   -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic linesearch
13946f263ca3SPeter Brune 
13956f263ca3SPeter Brune    Notes:
13961f8196a2SJed Brown    This is most relevant to the SNESLINESEARCHBASIC line search type since most line searches have a stopping criteria involving the norm.
13976f263ca3SPeter Brune 
13986f263ca3SPeter Brune    Level: intermediate
13996f263ca3SPeter Brune 
14001a4f838cSPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetNorms(), SNESLineSearchComputeNorms(), SNESLINESEARCHBASIC
14016f263ca3SPeter Brune @*/
14026f263ca3SPeter Brune PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg)
14036f263ca3SPeter Brune {
14046f263ca3SPeter Brune   PetscFunctionBegin;
14056f263ca3SPeter Brune   linesearch->norms = flg;
14066f263ca3SPeter Brune   PetscFunctionReturn(0);
14076f263ca3SPeter Brune }
14086f263ca3SPeter Brune 
1409f40b411bSPeter Brune /*@
1410f1c6b773SPeter Brune    SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context
1411f40b411bSPeter Brune 
1412f40b411bSPeter Brune    Input Parameters:
141378bcb3b5SPeter Brune .  linesearch - linesearch context
1414f40b411bSPeter Brune 
1415f40b411bSPeter Brune    Output Parameters:
14166232e825SPeter Brune +  X - Solution vector
14176232e825SPeter Brune .  F - Function vector
14186232e825SPeter Brune .  Y - Search direction vector
14196232e825SPeter Brune .  W - Solution work vector
14206232e825SPeter Brune -  G - Function work vector
14216232e825SPeter Brune 
14227bba9028SPeter Brune    Notes:
14237bba9028SPeter Brune    At the beginning of a line search application, X should contain a
14246232e825SPeter Brune    solution and the vector F the function computed at X.  At the end of the
14256232e825SPeter Brune    line search application, X should contain the new solution, and F the
14266232e825SPeter Brune    function evaluated at the new solution.
1427f40b411bSPeter Brune 
14282a7a6963SBarry Smith    These vectors are owned by the SNESLineSearch and should not be destroyed by the caller
14292a7a6963SBarry Smith 
143078bcb3b5SPeter Brune    Level: advanced
1431f40b411bSPeter Brune 
1432f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1433f40b411bSPeter Brune @*/
1434bf388a1fSBarry Smith PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G)
1435bf388a1fSBarry Smith {
14366a388c36SPeter Brune   PetscFunctionBegin;
1437f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
14386a388c36SPeter Brune   if (X) {
14396a388c36SPeter Brune     PetscValidPointer(X, 2);
14406a388c36SPeter Brune     *X = linesearch->vec_sol;
14416a388c36SPeter Brune   }
14426a388c36SPeter Brune   if (F) {
14436a388c36SPeter Brune     PetscValidPointer(F, 3);
14446a388c36SPeter Brune     *F = linesearch->vec_func;
14456a388c36SPeter Brune   }
14466a388c36SPeter Brune   if (Y) {
14476a388c36SPeter Brune     PetscValidPointer(Y, 4);
14486a388c36SPeter Brune     *Y = linesearch->vec_update;
14496a388c36SPeter Brune   }
14506a388c36SPeter Brune   if (W) {
14516a388c36SPeter Brune     PetscValidPointer(W, 5);
14526a388c36SPeter Brune     *W = linesearch->vec_sol_new;
14536a388c36SPeter Brune   }
14546a388c36SPeter Brune   if (G) {
14556a388c36SPeter Brune     PetscValidPointer(G, 6);
14566a388c36SPeter Brune     *G = linesearch->vec_func_new;
14576a388c36SPeter Brune   }
14586a388c36SPeter Brune   PetscFunctionReturn(0);
14596a388c36SPeter Brune }
14606a388c36SPeter Brune 
1461f40b411bSPeter Brune /*@
1462f1c6b773SPeter Brune    SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context
1463f40b411bSPeter Brune 
1464f40b411bSPeter Brune    Input Parameters:
146578bcb3b5SPeter Brune +  linesearch - linesearch context
14666232e825SPeter Brune .  X - Solution vector
14676232e825SPeter Brune .  F - Function vector
14686232e825SPeter Brune .  Y - Search direction vector
14696232e825SPeter Brune .  W - Solution work vector
14706232e825SPeter Brune -  G - Function work vector
1471f40b411bSPeter Brune 
147278bcb3b5SPeter Brune    Level: advanced
1473f40b411bSPeter Brune 
1474f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms(), SNESLineSearchGetVecs()
1475f40b411bSPeter Brune @*/
1476bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G)
1477bf388a1fSBarry Smith {
14786a388c36SPeter Brune   PetscFunctionBegin;
1479f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
14806a388c36SPeter Brune   if (X) {
14816a388c36SPeter Brune     PetscValidHeaderSpecific(X,VEC_CLASSID,2);
14826a388c36SPeter Brune     linesearch->vec_sol = X;
14836a388c36SPeter Brune   }
14846a388c36SPeter Brune   if (F) {
14856a388c36SPeter Brune     PetscValidHeaderSpecific(F,VEC_CLASSID,3);
14866a388c36SPeter Brune     linesearch->vec_func = F;
14876a388c36SPeter Brune   }
14886a388c36SPeter Brune   if (Y) {
14896a388c36SPeter Brune     PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
14906a388c36SPeter Brune     linesearch->vec_update = Y;
14916a388c36SPeter Brune   }
14926a388c36SPeter Brune   if (W) {
14936a388c36SPeter Brune     PetscValidHeaderSpecific(W,VEC_CLASSID,5);
14946a388c36SPeter Brune     linesearch->vec_sol_new = W;
14956a388c36SPeter Brune   }
14966a388c36SPeter Brune   if (G) {
14976a388c36SPeter Brune     PetscValidHeaderSpecific(G,VEC_CLASSID,6);
14986a388c36SPeter Brune     linesearch->vec_func_new = G;
14996a388c36SPeter Brune   }
15006a388c36SPeter Brune   PetscFunctionReturn(0);
15016a388c36SPeter Brune }
15026a388c36SPeter Brune 
1503e7058c64SPeter Brune /*@C
1504f1c6b773SPeter Brune    SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1505e7058c64SPeter Brune    SNES options in the database.
1506e7058c64SPeter Brune 
1507cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
1508e7058c64SPeter Brune 
1509e7058c64SPeter Brune    Input Parameters:
1510e7058c64SPeter Brune +  snes - the SNES context
1511e7058c64SPeter Brune -  prefix - the prefix to prepend to all option names
1512e7058c64SPeter Brune 
1513e7058c64SPeter Brune    Notes:
1514e7058c64SPeter Brune    A hyphen (-) must NOT be given at the beginning of the prefix name.
1515e7058c64SPeter Brune    The first character of all runtime options is AUTOMATICALLY the hyphen.
1516e7058c64SPeter Brune 
1517e7058c64SPeter Brune    Level: advanced
1518e7058c64SPeter Brune 
1519e7058c64SPeter Brune .seealso: SNESGetOptionsPrefix()
1520e7058c64SPeter Brune @*/
1521f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[])
1522e7058c64SPeter Brune {
1523e7058c64SPeter Brune   PetscErrorCode ierr;
1524e7058c64SPeter Brune 
1525e7058c64SPeter Brune   PetscFunctionBegin;
1526f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1527e7058c64SPeter Brune   ierr = PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1528e7058c64SPeter Brune   PetscFunctionReturn(0);
1529e7058c64SPeter Brune }
1530e7058c64SPeter Brune 
1531e7058c64SPeter Brune /*@C
1532f1c6b773SPeter Brune    SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all
1533f1c6b773SPeter Brune    SNESLineSearch options in the database.
1534e7058c64SPeter Brune 
1535e7058c64SPeter Brune    Not Collective
1536e7058c64SPeter Brune 
1537e7058c64SPeter Brune    Input Parameter:
1538cd7522eaSPeter Brune .  linesearch - the SNESLineSearch context
1539e7058c64SPeter Brune 
1540e7058c64SPeter Brune    Output Parameter:
1541e7058c64SPeter Brune .  prefix - pointer to the prefix string used
1542e7058c64SPeter Brune 
15438e557f58SPeter Brune    Notes:
15448e557f58SPeter Brune    On the fortran side, the user should pass in a string 'prefix' of
1545e7058c64SPeter Brune    sufficient length to hold the prefix.
1546e7058c64SPeter Brune 
1547e7058c64SPeter Brune    Level: advanced
1548e7058c64SPeter Brune 
1549e7058c64SPeter Brune .seealso: SNESAppendOptionsPrefix()
1550e7058c64SPeter Brune @*/
1551f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[])
1552e7058c64SPeter Brune {
1553e7058c64SPeter Brune   PetscErrorCode ierr;
1554e7058c64SPeter Brune 
1555e7058c64SPeter Brune   PetscFunctionBegin;
1556f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1557e7058c64SPeter Brune   ierr = PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1558e7058c64SPeter Brune   PetscFunctionReturn(0);
1559e7058c64SPeter Brune }
1560bf7f4e0aSPeter Brune 
15618d359177SBarry Smith /*@C
1562fa0ddf94SBarry Smith    SNESLineSearchSetWorkVecs - Gets work vectors for the line search.
1563f40b411bSPeter Brune 
1564f40b411bSPeter Brune    Input Parameter:
1565f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
1566f40b411bSPeter Brune -  nwork - the number of work vectors
1567f40b411bSPeter Brune 
1568f40b411bSPeter Brune    Level: developer
1569f40b411bSPeter Brune 
1570fa0ddf94SBarry Smith .seealso: SNESSetWorkVecs()
1571f40b411bSPeter Brune @*/
1572fa0ddf94SBarry Smith PetscErrorCode  SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork)
1573bf7f4e0aSPeter Brune {
1574bf7f4e0aSPeter Brune   PetscErrorCode ierr;
1575bf388a1fSBarry Smith 
1576bf7f4e0aSPeter Brune   PetscFunctionBegin;
1577bf7f4e0aSPeter Brune   if (linesearch->vec_sol) {
1578bf7f4e0aSPeter Brune     ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr);
15798d359177SBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
1580bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1581bf7f4e0aSPeter Brune }
1582bf7f4e0aSPeter Brune 
1583f40b411bSPeter Brune /*@
1584422a814eSBarry Smith    SNESLineSearchGetReason - Gets the success/failure status of the last line search application
1585f40b411bSPeter Brune 
1586f40b411bSPeter Brune    Input Parameters:
158778bcb3b5SPeter Brune .  linesearch - linesearch context
1588f40b411bSPeter Brune 
1589f40b411bSPeter Brune    Output Parameters:
1590422a814eSBarry Smith .  result - The success or failure status
1591f40b411bSPeter Brune 
1592cd7522eaSPeter Brune    Notes:
1593c98378a5SBarry Smith    This is typically called after SNESLineSearchApply() in order to determine if the line-search failed
1594cd7522eaSPeter Brune    (and set the SNES convergence accordingly).
1595cd7522eaSPeter Brune 
1596f40b411bSPeter Brune    Level: intermediate
1597f40b411bSPeter Brune 
1598422a814eSBarry Smith .seealso: SNESLineSearchSetReason(), SNESLineSearchReason
1599f40b411bSPeter Brune @*/
1600422a814eSBarry Smith PetscErrorCode  SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result)
1601bf7f4e0aSPeter Brune {
1602bf7f4e0aSPeter Brune   PetscFunctionBegin;
1603f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1604422a814eSBarry Smith   PetscValidPointer(result, 2);
1605422a814eSBarry Smith   *result = linesearch->result;
1606bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1607bf7f4e0aSPeter Brune }
1608bf7f4e0aSPeter Brune 
1609f40b411bSPeter Brune /*@
1610422a814eSBarry Smith    SNESLineSearchSetReason - Sets the success/failure status of the last line search application
1611f40b411bSPeter Brune 
1612f40b411bSPeter Brune    Input Parameters:
161378bcb3b5SPeter Brune +  linesearch - linesearch context
1614422a814eSBarry Smith -  result - The success or failure status
1615f40b411bSPeter Brune 
1616cd7522eaSPeter Brune    Notes:
1617cd7522eaSPeter Brune    This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set
1618cd7522eaSPeter Brune    the success or failure of the line search method.
1619cd7522eaSPeter Brune 
162078bcb3b5SPeter Brune    Level: developer
1621f40b411bSPeter Brune 
1622422a814eSBarry Smith .seealso: SNESLineSearchGetSResult()
1623f40b411bSPeter Brune @*/
1624422a814eSBarry Smith PetscErrorCode  SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result)
16256a388c36SPeter Brune {
16266a388c36SPeter Brune   PetscFunctionBegin;
16275fd66863SKarl Rupp   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1628422a814eSBarry Smith   linesearch->result = result;
16296a388c36SPeter Brune   PetscFunctionReturn(0);
16306a388c36SPeter Brune }
16316a388c36SPeter Brune 
16329bd66eb0SPeter Brune /*@C
1633f1c6b773SPeter Brune    SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.
16349bd66eb0SPeter Brune 
16359bd66eb0SPeter Brune    Input Parameters:
16369bd66eb0SPeter Brune +  snes - nonlinear context obtained from SNESCreate()
16379bd66eb0SPeter Brune .  projectfunc - function for projecting the function to the bounds
16389bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
16399bd66eb0SPeter Brune 
16409bd66eb0SPeter Brune    Logically Collective on SNES
16419bd66eb0SPeter Brune 
16429bd66eb0SPeter Brune    Calling sequence of projectfunc:
16439bd66eb0SPeter Brune .vb
16449bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X)
16459bd66eb0SPeter Brune .ve
16469bd66eb0SPeter Brune 
16479bd66eb0SPeter Brune     Input parameters for projectfunc:
16489bd66eb0SPeter Brune +   snes - nonlinear context
16499bd66eb0SPeter Brune -   X - current solution
16509bd66eb0SPeter Brune 
1651cd7522eaSPeter Brune     Output parameters for projectfunc:
16529bd66eb0SPeter Brune .   X - Projected solution
16539bd66eb0SPeter Brune 
16549bd66eb0SPeter Brune    Calling sequence of normfunc:
16559bd66eb0SPeter Brune .vb
16569bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm)
16579bd66eb0SPeter Brune .ve
16589bd66eb0SPeter Brune 
1659cd7522eaSPeter Brune     Input parameters for normfunc:
16609bd66eb0SPeter Brune +   snes - nonlinear context
16619bd66eb0SPeter Brune .   X - current solution
16629bd66eb0SPeter Brune -   F - current residual
16639bd66eb0SPeter Brune 
1664cd7522eaSPeter Brune     Output parameters for normfunc:
16659bd66eb0SPeter Brune .   fnorm - VI-specific norm of the function
16669bd66eb0SPeter Brune 
1667cd7522eaSPeter Brune     Notes:
1668cd7522eaSPeter Brune     The VI solvers require projection of the solution to the feasible set.  projectfunc should implement this.
1669cd7522eaSPeter Brune 
1670cd7522eaSPeter Brune     The VI solvers require special evaluation of the function norm such that the norm is only calculated
1671cd7522eaSPeter Brune     on the inactive set.  This should be implemented by normfunc.
16729bd66eb0SPeter Brune 
16739bd66eb0SPeter Brune     Level: developer
16749bd66eb0SPeter Brune 
1675f1c6b773SPeter Brune .seealso: SNESLineSearchGetVIFunctions(), SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck()
16769bd66eb0SPeter Brune @*/
167725acbd8eSLisandro Dalcin PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc)
16789bd66eb0SPeter Brune {
16799bd66eb0SPeter Brune   PetscFunctionBegin;
1680f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
16819bd66eb0SPeter Brune   if (projectfunc) linesearch->ops->viproject = projectfunc;
16829bd66eb0SPeter Brune   if (normfunc) linesearch->ops->vinorm = normfunc;
16839bd66eb0SPeter Brune   PetscFunctionReturn(0);
16849bd66eb0SPeter Brune }
16859bd66eb0SPeter Brune 
16869bd66eb0SPeter Brune /*@C
1687f1c6b773SPeter Brune    SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.
16889bd66eb0SPeter Brune 
16899bd66eb0SPeter Brune    Input Parameters:
1690907376e6SBarry Smith .  linesearch - the line search context, obtain with SNESGetLineSearch()
16919bd66eb0SPeter Brune 
16929bd66eb0SPeter Brune    Output Parameters:
16939bd66eb0SPeter Brune +  projectfunc - function for projecting the function to the bounds
16949bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
16959bd66eb0SPeter Brune 
16969bd66eb0SPeter Brune    Logically Collective on SNES
16979bd66eb0SPeter Brune 
16989bd66eb0SPeter Brune     Level: developer
16999bd66eb0SPeter Brune 
1700f1c6b773SPeter Brune .seealso: SNESLineSearchSetVIFunctions(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
17019bd66eb0SPeter Brune @*/
170225acbd8eSLisandro Dalcin PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc)
17039bd66eb0SPeter Brune {
17049bd66eb0SPeter Brune   PetscFunctionBegin;
17059bd66eb0SPeter Brune   if (projectfunc) *projectfunc = linesearch->ops->viproject;
17069bd66eb0SPeter Brune   if (normfunc) *normfunc = linesearch->ops->vinorm;
17079bd66eb0SPeter Brune   PetscFunctionReturn(0);
17089bd66eb0SPeter Brune }
17099bd66eb0SPeter Brune 
1710bf7f4e0aSPeter Brune /*@C
17111c84c290SBarry Smith   SNESLineSearchRegister - See SNESLineSearchRegister()
1712bf7f4e0aSPeter Brune 
1713bf7f4e0aSPeter Brune   Level: advanced
1714bf7f4e0aSPeter Brune @*/
1715bdf89e91SBarry Smith PetscErrorCode  SNESLineSearchRegister(const char sname[],PetscErrorCode (*function)(SNESLineSearch))
1716bf7f4e0aSPeter Brune {
1717bf7f4e0aSPeter Brune   PetscErrorCode ierr;
1718bf7f4e0aSPeter Brune 
1719bf7f4e0aSPeter Brune   PetscFunctionBegin;
17201d36bdfdSBarry Smith   ierr = SNESInitializePackage();CHKERRQ(ierr);
1721a240a19fSJed Brown   ierr = PetscFunctionListAdd(&SNESLineSearchList,sname,function);CHKERRQ(ierr);
1722bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1723bf7f4e0aSPeter Brune }
1724