xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision d8d19677bbccf95218448bee62e6b87f4513e133)
1af0996ceSBarry Smith #include <petsc/private/linesearchimpl.h> /*I "petscsnes.h" I*/
2bf7f4e0aSPeter Brune 
3f1c6b773SPeter Brune PetscBool         SNESLineSearchRegisterAllCalled = PETSC_FALSE;
40298fd71SBarry Smith PetscFunctionList SNESLineSearchList              = NULL;
5bf7f4e0aSPeter Brune 
6f1c6b773SPeter Brune PetscClassId  SNESLINESEARCH_CLASSID;
757a83faaSBarry Smith PetscLogEvent SNESLINESEARCH_Apply;
8bf7f4e0aSPeter Brune 
9dcf2fd19SBarry Smith /*@
10dcf2fd19SBarry Smith    SNESLineSearchMonitorCancel - Clears all the monitor functions for a SNESLineSearch object.
11dcf2fd19SBarry Smith 
12dcf2fd19SBarry Smith    Logically Collective on SNESLineSearch
13dcf2fd19SBarry Smith 
14dcf2fd19SBarry Smith    Input Parameters:
15dcf2fd19SBarry Smith .  ls - the SNESLineSearch context
16dcf2fd19SBarry Smith 
17dcf2fd19SBarry Smith    Options Database Key:
18dcf2fd19SBarry Smith .  -snes_linesearch_monitor_cancel - cancels all monitors that have been hardwired
19dcf2fd19SBarry Smith     into a code by calls to SNESLineSearchMonitorSet(), but does not cancel those
20dcf2fd19SBarry Smith     set via the options database
21dcf2fd19SBarry Smith 
22dcf2fd19SBarry Smith    Notes:
23dcf2fd19SBarry Smith    There is no way to clear one specific monitor from a SNESLineSearch object.
24dcf2fd19SBarry Smith 
25dcf2fd19SBarry Smith    This does not clear the monitor set with SNESLineSearchSetDefaultMonitor() use SNESLineSearchSetDefaultMonitor(ls,NULL) to cancel
26dcf2fd19SBarry Smith    that one.
27dcf2fd19SBarry Smith 
28dcf2fd19SBarry Smith    Level: intermediate
29dcf2fd19SBarry Smith 
3084238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchMonitorDefault(), SNESLineSearchMonitorSet()
31dcf2fd19SBarry Smith @*/
32dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchMonitorCancel(SNESLineSearch ls)
33dcf2fd19SBarry Smith {
34dcf2fd19SBarry Smith   PetscErrorCode ierr;
35dcf2fd19SBarry Smith   PetscInt       i;
36dcf2fd19SBarry Smith 
37dcf2fd19SBarry Smith   PetscFunctionBegin;
38dcf2fd19SBarry Smith   PetscValidHeaderSpecific(ls,SNESLINESEARCH_CLASSID,1);
39dcf2fd19SBarry Smith   for (i=0; i<ls->numbermonitors; i++) {
40dcf2fd19SBarry Smith     if (ls->monitordestroy[i]) {
41dcf2fd19SBarry Smith       ierr = (*ls->monitordestroy[i])(&ls->monitorcontext[i]);CHKERRQ(ierr);
42dcf2fd19SBarry Smith     }
43dcf2fd19SBarry Smith   }
44dcf2fd19SBarry Smith   ls->numbermonitors = 0;
45dcf2fd19SBarry Smith   PetscFunctionReturn(0);
46dcf2fd19SBarry Smith }
47dcf2fd19SBarry Smith 
48dcf2fd19SBarry Smith /*@
49dcf2fd19SBarry Smith    SNESLineSearchMonitor - runs the user provided monitor routines, if they exist
50dcf2fd19SBarry Smith 
51dcf2fd19SBarry Smith    Collective on SNES
52dcf2fd19SBarry Smith 
53dcf2fd19SBarry Smith    Input Parameters:
54dcf2fd19SBarry Smith .  ls - the linesearch object
55dcf2fd19SBarry Smith 
56dcf2fd19SBarry Smith    Notes:
57dcf2fd19SBarry Smith    This routine is called by the SNES implementations.
58dcf2fd19SBarry Smith    It does not typically need to be called by the user.
59dcf2fd19SBarry Smith 
60dcf2fd19SBarry Smith    Level: developer
61dcf2fd19SBarry Smith 
6284238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchMonitorSet()
63dcf2fd19SBarry Smith @*/
64dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchMonitor(SNESLineSearch ls)
65dcf2fd19SBarry Smith {
66dcf2fd19SBarry Smith   PetscErrorCode ierr;
67dcf2fd19SBarry Smith   PetscInt       i,n = ls->numbermonitors;
68dcf2fd19SBarry Smith 
69dcf2fd19SBarry Smith   PetscFunctionBegin;
70dcf2fd19SBarry Smith   for (i=0; i<n; i++) {
71dcf2fd19SBarry Smith     ierr = (*ls->monitorftns[i])(ls,ls->monitorcontext[i]);CHKERRQ(ierr);
72dcf2fd19SBarry Smith   }
73dcf2fd19SBarry Smith   PetscFunctionReturn(0);
74dcf2fd19SBarry Smith }
75dcf2fd19SBarry Smith 
76dcf2fd19SBarry Smith /*@C
77dcf2fd19SBarry Smith    SNESLineSearchMonitorSet - Sets an ADDITIONAL function that is to be used at every
78dcf2fd19SBarry Smith    iteration of the nonlinear solver to display the iteration's
79dcf2fd19SBarry Smith    progress.
80dcf2fd19SBarry Smith 
81dcf2fd19SBarry Smith    Logically Collective on SNESLineSearch
82dcf2fd19SBarry Smith 
83dcf2fd19SBarry Smith    Input Parameters:
84dcf2fd19SBarry Smith +  ls - the SNESLineSearch context
85dcf2fd19SBarry Smith .  f - the monitor function
86dcf2fd19SBarry Smith .  mctx - [optional] user-defined context for private data for the
87dcf2fd19SBarry Smith           monitor routine (use NULL if no context is desired)
88dcf2fd19SBarry Smith -  monitordestroy - [optional] routine that frees monitor context
89dcf2fd19SBarry Smith           (may be NULL)
90dcf2fd19SBarry Smith 
91dcf2fd19SBarry Smith    Notes:
92dcf2fd19SBarry Smith    Several different monitoring routines may be set by calling
93dcf2fd19SBarry Smith    SNESLineSearchMonitorSet() multiple times; all will be called in the
94dcf2fd19SBarry Smith    order in which they were set.
95dcf2fd19SBarry Smith 
9695452b02SPatrick Sanan    Fortran Notes:
9795452b02SPatrick Sanan     Only a single monitor function can be set for each SNESLineSearch object
98dcf2fd19SBarry Smith 
99dcf2fd19SBarry Smith    Level: intermediate
100dcf2fd19SBarry Smith 
10184238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchMonitorDefault(), SNESLineSearchMonitorCancel()
102dcf2fd19SBarry Smith @*/
103dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchMonitorSet(SNESLineSearch ls,PetscErrorCode (*f)(SNESLineSearch,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
104dcf2fd19SBarry Smith {
10578064530SBarry Smith   PetscErrorCode ierr;
10678064530SBarry Smith   PetscInt       i;
10778064530SBarry Smith   PetscBool      identical;
10878064530SBarry Smith 
109dcf2fd19SBarry Smith   PetscFunctionBegin;
110dcf2fd19SBarry Smith   PetscValidHeaderSpecific(ls,SNESLINESEARCH_CLASSID,1);
11178064530SBarry Smith   for (i=0; i<ls->numbermonitors;i++) {
11278064530SBarry Smith     ierr = PetscMonitorCompare((PetscErrorCode (*)(void))f,mctx,monitordestroy,(PetscErrorCode (*)(void))ls->monitorftns[i],ls->monitorcontext[i],ls->monitordestroy[i],&identical);CHKERRQ(ierr);
11378064530SBarry Smith     if (identical) PetscFunctionReturn(0);
11478064530SBarry Smith   }
115dcf2fd19SBarry Smith   if (ls->numbermonitors >= MAXSNESLSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
116dcf2fd19SBarry Smith   ls->monitorftns[ls->numbermonitors]          = f;
117dcf2fd19SBarry Smith   ls->monitordestroy[ls->numbermonitors]   = monitordestroy;
118dcf2fd19SBarry Smith   ls->monitorcontext[ls->numbermonitors++] = (void*)mctx;
119dcf2fd19SBarry Smith   PetscFunctionReturn(0);
120dcf2fd19SBarry Smith }
121dcf2fd19SBarry Smith 
122dcf2fd19SBarry Smith /*@C
123dcf2fd19SBarry Smith    SNESLineSearchMonitorSolutionUpdate - Monitors each update a new function value the linesearch tries
124dcf2fd19SBarry Smith 
125dcf2fd19SBarry Smith    Collective on SNESLineSearch
126dcf2fd19SBarry Smith 
127dcf2fd19SBarry Smith    Input Parameters:
128dcf2fd19SBarry Smith +  ls - the SNES linesearch object
129d12e167eSBarry Smith -  vf - the context for the monitor, in this case it is an ASCII PetscViewer and format
130dcf2fd19SBarry Smith 
131dcf2fd19SBarry Smith    Level: intermediate
132dcf2fd19SBarry Smith 
13384238204SBarry Smith .seealso: SNESGetLineSearch(), SNESMonitorSet(), SNESMonitorSolution()
134dcf2fd19SBarry Smith @*/
135d12e167eSBarry Smith PetscErrorCode  SNESLineSearchMonitorSolutionUpdate(SNESLineSearch ls,PetscViewerAndFormat *vf)
136dcf2fd19SBarry Smith {
137dcf2fd19SBarry Smith   PetscErrorCode ierr;
138d12e167eSBarry Smith   PetscViewer    viewer = vf->viewer;
139dcf2fd19SBarry Smith   Vec            Y,W,G;
140dcf2fd19SBarry Smith 
141dcf2fd19SBarry Smith   PetscFunctionBegin;
142dcf2fd19SBarry Smith   ierr = SNESLineSearchGetVecs(ls,NULL,NULL,&Y,&W,&G);CHKERRQ(ierr);
143d12e167eSBarry Smith   ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr);
144dcf2fd19SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"LineSearch attempted update to solution \n");CHKERRQ(ierr);
145dcf2fd19SBarry Smith   ierr = VecView(Y,viewer);CHKERRQ(ierr);
146dcf2fd19SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"LineSearch attempted new solution \n");CHKERRQ(ierr);
147dcf2fd19SBarry Smith   ierr = VecView(W,viewer);CHKERRQ(ierr);
148dcf2fd19SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"LineSearch attempted updated function value\n");CHKERRQ(ierr);
149dcf2fd19SBarry Smith   ierr = VecView(G,viewer);CHKERRQ(ierr);
150d12e167eSBarry Smith   ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
151dcf2fd19SBarry Smith   PetscFunctionReturn(0);
152dcf2fd19SBarry Smith }
153dcf2fd19SBarry Smith 
154f40b411bSPeter Brune /*@
155cd7522eaSPeter Brune    SNESLineSearchCreate - Creates the line search context.
156f40b411bSPeter Brune 
157cd7522eaSPeter Brune    Logically Collective on Comm
158f40b411bSPeter Brune 
159f40b411bSPeter Brune    Input Parameters:
160cd7522eaSPeter Brune .  comm - MPI communicator for the line search (typically from the associated SNES context).
161f40b411bSPeter Brune 
162f40b411bSPeter Brune    Output Parameters:
1638e557f58SPeter Brune .  outlinesearch - the new linesearch context
164f40b411bSPeter Brune 
165162e0bf5SPeter Brune    Level: developer
166162e0bf5SPeter Brune 
167162e0bf5SPeter Brune    Notes:
168162e0bf5SPeter Brune    The preferred calling sequence for users is to use SNESGetLineSearch() to acquire the SNESLineSearch instance
169162e0bf5SPeter Brune    already associated with the SNES.  This function is for developer use.
170f40b411bSPeter Brune 
171162e0bf5SPeter Brune .seealso: LineSearchDestroy(), SNESGetLineSearch()
172f40b411bSPeter Brune @*/
173f40b411bSPeter Brune 
174bf388a1fSBarry Smith PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch)
175bf388a1fSBarry Smith {
176bf7f4e0aSPeter Brune   PetscErrorCode ierr;
177f1c6b773SPeter Brune   SNESLineSearch linesearch;
178bf388a1fSBarry Smith 
179bf7f4e0aSPeter Brune   PetscFunctionBegin;
180ea5d4fccSPeter Brune   PetscValidPointer(outlinesearch,2);
181f34a81feSMatthew G. Knepley   ierr = SNESInitializePackage();CHKERRQ(ierr);
1820298fd71SBarry Smith   *outlinesearch = NULL;
183f5af7f23SKarl Rupp 
18473107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(linesearch,SNESLINESEARCH_CLASSID, "SNESLineSearch","Linesearch","SNESLineSearch",comm,SNESLineSearchDestroy,SNESLineSearchView);CHKERRQ(ierr);
185bf7f4e0aSPeter Brune 
1860298fd71SBarry Smith   linesearch->vec_sol_new  = NULL;
1870298fd71SBarry Smith   linesearch->vec_func_new = NULL;
1880298fd71SBarry Smith   linesearch->vec_sol      = NULL;
1890298fd71SBarry Smith   linesearch->vec_func     = NULL;
1900298fd71SBarry Smith   linesearch->vec_update   = NULL;
1919bd66eb0SPeter Brune 
192bf7f4e0aSPeter Brune   linesearch->lambda       = 1.0;
193bf7f4e0aSPeter Brune   linesearch->fnorm        = 1.0;
194bf7f4e0aSPeter Brune   linesearch->ynorm        = 1.0;
195bf7f4e0aSPeter Brune   linesearch->xnorm        = 1.0;
196422a814eSBarry Smith   linesearch->result       = SNES_LINESEARCH_SUCCEEDED;
197bf7f4e0aSPeter Brune   linesearch->norms        = PETSC_TRUE;
198bf7f4e0aSPeter Brune   linesearch->keeplambda   = PETSC_FALSE;
199bf7f4e0aSPeter Brune   linesearch->damping      = 1.0;
200bf7f4e0aSPeter Brune   linesearch->maxstep      = 1e8;
201bf7f4e0aSPeter Brune   linesearch->steptol      = 1e-12;
202516fe3c3SPeter Brune   linesearch->rtol         = 1e-8;
203516fe3c3SPeter Brune   linesearch->atol         = 1e-15;
204516fe3c3SPeter Brune   linesearch->ltol         = 1e-8;
2050298fd71SBarry Smith   linesearch->precheckctx  = NULL;
2060298fd71SBarry Smith   linesearch->postcheckctx = NULL;
20759405d5eSPeter Brune   linesearch->max_its      = 1;
208bf7f4e0aSPeter Brune   linesearch->setupcalled  = PETSC_FALSE;
2093add74b1SBarry Smith   linesearch->monitor      = NULL;
210bf7f4e0aSPeter Brune   *outlinesearch           = linesearch;
211bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
212bf7f4e0aSPeter Brune }
213bf7f4e0aSPeter Brune 
214f40b411bSPeter Brune /*@
21578bcb3b5SPeter Brune    SNESLineSearchSetUp - Prepares the line search for being applied by allocating
21678bcb3b5SPeter Brune    any required vectors.
217f40b411bSPeter Brune 
218cd7522eaSPeter Brune    Collective on SNESLineSearch
219f40b411bSPeter Brune 
220f40b411bSPeter Brune    Input Parameters:
221f40b411bSPeter Brune .  linesearch - The LineSearch instance.
222f40b411bSPeter Brune 
223cd7522eaSPeter Brune    Notes:
224f190f2fcSBarry Smith    For most cases, this needn't be called by users or outside of SNESLineSearchApply().
225cd7522eaSPeter Brune    The only current case where this is called outside of this is for the VI
22678bcb3b5SPeter Brune    solvers, which modify the solution and work vectors before the first call
227cd7522eaSPeter Brune    of SNESLineSearchApply, requiring the SNESLineSearch work vectors to be
228cd7522eaSPeter Brune    allocated upfront.
229cd7522eaSPeter Brune 
23078bcb3b5SPeter Brune    Level: advanced
231f40b411bSPeter Brune 
23284238204SBarry Smith .seealso: SNESGetLineSearch(), SNESLineSearchReset()
233f40b411bSPeter Brune @*/
234f40b411bSPeter Brune 
235bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch)
236bf388a1fSBarry Smith {
237bf7f4e0aSPeter Brune   PetscErrorCode ierr;
238bf388a1fSBarry Smith 
239bf7f4e0aSPeter Brune   PetscFunctionBegin;
240bf7f4e0aSPeter Brune   if (!((PetscObject)linesearch)->type_name) {
2411a4f838cSPeter Brune     ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr);
242bf7f4e0aSPeter Brune   }
243bf7f4e0aSPeter Brune   if (!linesearch->setupcalled) {
2449bd66eb0SPeter Brune     if (!linesearch->vec_sol_new) {
245bf7f4e0aSPeter Brune       ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new);CHKERRQ(ierr);
2469bd66eb0SPeter Brune     }
2479bd66eb0SPeter Brune     if (!linesearch->vec_func_new) {
2489bd66eb0SPeter Brune       ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new);CHKERRQ(ierr);
2499bd66eb0SPeter Brune     }
250bf7f4e0aSPeter Brune     if (linesearch->ops->setup) {
251bf7f4e0aSPeter Brune       ierr = (*linesearch->ops->setup)(linesearch);CHKERRQ(ierr);
252bf7f4e0aSPeter Brune     }
253ed07d7d7SPeter Brune     if (!linesearch->ops->snesfunc) {ierr = SNESLineSearchSetFunction(linesearch,SNESComputeFunction);CHKERRQ(ierr);}
254bf7f4e0aSPeter Brune     linesearch->lambda      = linesearch->damping;
255bf7f4e0aSPeter Brune     linesearch->setupcalled = PETSC_TRUE;
256bf7f4e0aSPeter Brune   }
257bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
258bf7f4e0aSPeter Brune }
259bf7f4e0aSPeter Brune 
260f40b411bSPeter Brune /*@
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);
602064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(Y,VEC_CLASSID,5);
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);
6489e5d0892SLisandro Dalcin   if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = NULL; 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   *monitor = linesearch->monitor;
7126a388c36SPeter Brune   PetscFunctionReturn(0);
7136a388c36SPeter Brune }
7146a388c36SPeter Brune 
715dcf2fd19SBarry Smith /*@C
716dcf2fd19SBarry Smith    SNESLineSearchMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
717dcf2fd19SBarry Smith 
718dcf2fd19SBarry Smith    Collective on SNESLineSearch
719dcf2fd19SBarry Smith 
720dcf2fd19SBarry Smith    Input Parameters:
721dcf2fd19SBarry Smith +  ls - LineSearch object you wish to monitor
722dcf2fd19SBarry Smith .  name - the monitor type one is seeking
723dcf2fd19SBarry Smith .  help - message indicating what monitoring is done
724dcf2fd19SBarry Smith .  manual - manual page for the monitor
725dcf2fd19SBarry Smith .  monitor - the monitor function
726dcf2fd19SBarry 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
727dcf2fd19SBarry Smith 
728dcf2fd19SBarry Smith    Level: developer
729dcf2fd19SBarry Smith 
730dcf2fd19SBarry Smith .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
731dcf2fd19SBarry Smith           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
732dcf2fd19SBarry Smith           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
733dcf2fd19SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
734dcf2fd19SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
735dcf2fd19SBarry Smith           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
736dcf2fd19SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
737dcf2fd19SBarry Smith @*/
738d12e167eSBarry Smith PetscErrorCode  SNESLineSearchMonitorSetFromOptions(SNESLineSearch ls,const char name[],const char help[], const char manual[],PetscErrorCode (*monitor)(SNESLineSearch,PetscViewerAndFormat*),PetscErrorCode (*monitorsetup)(SNESLineSearch,PetscViewerAndFormat*))
739dcf2fd19SBarry Smith {
740dcf2fd19SBarry Smith   PetscErrorCode    ierr;
741dcf2fd19SBarry Smith   PetscViewer       viewer;
742dcf2fd19SBarry Smith   PetscViewerFormat format;
743dcf2fd19SBarry Smith   PetscBool         flg;
744dcf2fd19SBarry Smith 
745dcf2fd19SBarry Smith   PetscFunctionBegin;
74616413a6aSBarry Smith   ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)ls),((PetscObject) ls)->options,((PetscObject)ls)->prefix,name,&viewer,&format,&flg);CHKERRQ(ierr);
747dcf2fd19SBarry Smith   if (flg) {
748d12e167eSBarry Smith     PetscViewerAndFormat *vf;
749d12e167eSBarry Smith     ierr = PetscViewerAndFormatCreate(viewer,format,&vf);CHKERRQ(ierr);
750d12e167eSBarry Smith     ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr);
751dcf2fd19SBarry Smith     if (monitorsetup) {
752d12e167eSBarry Smith       ierr = (*monitorsetup)(ls,vf);CHKERRQ(ierr);
753dcf2fd19SBarry Smith     }
754d12e167eSBarry Smith     ierr = SNESLineSearchMonitorSet(ls,(PetscErrorCode (*)(SNESLineSearch,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);CHKERRQ(ierr);
755dcf2fd19SBarry Smith   }
756dcf2fd19SBarry Smith   PetscFunctionReturn(0);
757dcf2fd19SBarry Smith }
758dcf2fd19SBarry Smith 
759f40b411bSPeter Brune /*@
760f1c6b773SPeter Brune    SNESLineSearchSetFromOptions - Sets options for the line search
761f40b411bSPeter Brune 
762f40b411bSPeter Brune    Input Parameters:
7638e557f58SPeter Brune .  linesearch - linesearch context
764f40b411bSPeter Brune 
765f40b411bSPeter Brune    Options Database Keys:
766d4c6564cSPatrick Farrell + -snes_linesearch_type <type> - basic, bt, l2, cp, nleqerr, shell
7673c7d6663SPeter Brune . -snes_linesearch_order <order> - 1, 2, 3.  Most types only support certain orders (bt supports 2 or 3)
7681fe24845SBarry Smith . -snes_linesearch_norms   - Turn on/off the linesearch norms for the basic linesearch typem (SNESLineSearchSetComputeNorms())
76971eef1aeSPeter Brune . -snes_linesearch_minlambda - The minimum step length
7701a9b3a06SPeter Brune . -snes_linesearch_maxstep - The maximum step size
7711a9b3a06SPeter Brune . -snes_linesearch_rtol - Relative tolerance for iterative line searches
7721a9b3a06SPeter Brune . -snes_linesearch_atol - Absolute tolerance for iterative line searches
7731a9b3a06SPeter Brune . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches
7741a9b3a06SPeter Brune . -snes_linesearch_max_it - The number of iterations for iterative line searches
775dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches
776dcf2fd19SBarry Smith . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine
7778e557f58SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter
778cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess.
7791a9b3a06SPeter Brune . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method
780d8d34be6SBarry Smith - -snes_linesearch_precheck_picard_angle - Angle used in Picard precheck method
781f40b411bSPeter Brune 
782f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
783f40b411bSPeter Brune 
784f40b411bSPeter Brune    Level: intermediate
785f40b411bSPeter Brune 
7861fe24845SBarry Smith .seealso: SNESLineSearchCreate(), SNESLineSearchSetOrder(), SNESLineSearchSetType(), SNESLineSearchSetTolerances(), SNESLineSearchSetDamping(), SNESLineSearchPreCheckPicard(),
7871fe24845SBarry Smith           SNESLineSearchType, SNESLineSearchSetComputeNorms()
788f40b411bSPeter Brune @*/
789bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch)
790bf388a1fSBarry Smith {
791bf7f4e0aSPeter Brune   PetscErrorCode    ierr;
7921a4f838cSPeter Brune   const char        *deft = SNESLINESEARCHBASIC;
793bf7f4e0aSPeter Brune   char              type[256];
794bf7f4e0aSPeter Brune   PetscBool         flg, set;
795dcf2fd19SBarry Smith   PetscViewer       viewer;
796bf388a1fSBarry Smith 
797bf7f4e0aSPeter Brune   PetscFunctionBegin;
7980f51fdf8SToby Isaac   ierr = SNESLineSearchRegisterAll();CHKERRQ(ierr);
799bf7f4e0aSPeter Brune 
800bf7f4e0aSPeter Brune   ierr = PetscObjectOptionsBegin((PetscObject)linesearch);CHKERRQ(ierr);
801f5af7f23SKarl Rupp   if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name;
802a264d7a6SBarry Smith   ierr = PetscOptionsFList("-snes_linesearch_type","Linesearch type","SNESLineSearchSetType",SNESLineSearchList,deft,type,256,&flg);CHKERRQ(ierr);
803bf7f4e0aSPeter Brune   if (flg) {
804f1c6b773SPeter Brune     ierr = SNESLineSearchSetType(linesearch,type);CHKERRQ(ierr);
805bf7f4e0aSPeter Brune   } else if (!((PetscObject)linesearch)->type_name) {
806f1c6b773SPeter Brune     ierr = SNESLineSearchSetType(linesearch,deft);CHKERRQ(ierr);
807bf7f4e0aSPeter Brune   }
808bf7f4e0aSPeter Brune 
80916413a6aSBarry Smith   ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)linesearch),((PetscObject) linesearch)->options,((PetscObject)linesearch)->prefix,"-snes_linesearch_monitor",&viewer,NULL,&set);CHKERRQ(ierr);
810dcf2fd19SBarry Smith   if (set) {
811dcf2fd19SBarry Smith     ierr = SNESLineSearchSetDefaultMonitor(linesearch,viewer);CHKERRQ(ierr);
812dcf2fd19SBarry Smith     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
813dcf2fd19SBarry Smith   }
814dcf2fd19SBarry Smith   ierr = SNESLineSearchMonitorSetFromOptions(linesearch,"-snes_linesearch_monitor_solution_update","View correction at each iteration","SNESLineSearchMonitorSolutionUpdate",SNESLineSearchMonitorSolutionUpdate,NULL);CHKERRQ(ierr);
815bf7f4e0aSPeter Brune 
8161a9b3a06SPeter Brune   /* tolerances */
81794ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_minlambda","Minimum step length","SNESLineSearchSetTolerances",linesearch->steptol,&linesearch->steptol,NULL);CHKERRQ(ierr);
81894ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_maxstep","Maximum step size","SNESLineSearchSetTolerances",linesearch->maxstep,&linesearch->maxstep,NULL);CHKERRQ(ierr);
81994ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_rtol","Relative tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->rtol,&linesearch->rtol,NULL);CHKERRQ(ierr);
82094ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_atol","Absolute tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->atol,&linesearch->atol,NULL);CHKERRQ(ierr);
82194ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_ltol","Change in lambda tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->ltol,&linesearch->ltol,NULL);CHKERRQ(ierr);
82294ae4db5SBarry Smith   ierr = PetscOptionsInt("-snes_linesearch_max_it","Maximum iterations for iterative line searches","SNESLineSearchSetTolerances",linesearch->max_its,&linesearch->max_its,NULL);CHKERRQ(ierr);
8237a35526eSPeter Brune 
8241a9b3a06SPeter Brune   /* damping parameters */
82594ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_damping","Line search damping and initial step guess","SNESLineSearchSetDamping",linesearch->damping,&linesearch->damping,NULL);CHKERRQ(ierr);
8261a9b3a06SPeter Brune 
82794ae4db5SBarry Smith   ierr = PetscOptionsBool("-snes_linesearch_keeplambda","Use previous lambda as damping","SNESLineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,NULL);CHKERRQ(ierr);
8281a9b3a06SPeter Brune 
8291a9b3a06SPeter Brune   /* precheck */
8307a35526eSPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_precheck_picard","Use a correction that sometimes improves convergence of Picard iteration","SNESLineSearchPreCheckPicard",flg,&flg,&set);CHKERRQ(ierr);
8317a35526eSPeter Brune   if (set) {
8327a35526eSPeter Brune     if (flg) {
8337a35526eSPeter Brune       linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */
834f5af7f23SKarl Rupp 
8357a35526eSPeter Brune       ierr = PetscOptionsReal("-snes_linesearch_precheck_picard_angle","Maximum angle at which to activate the correction",
8360298fd71SBarry Smith                               "none",linesearch->precheck_picard_angle,&linesearch->precheck_picard_angle,NULL);CHKERRQ(ierr);
8377a35526eSPeter Brune       ierr = SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&linesearch->precheck_picard_angle);CHKERRQ(ierr);
8387a35526eSPeter Brune     } else {
8390298fd71SBarry Smith       ierr = SNESLineSearchSetPreCheck(linesearch,NULL,NULL);CHKERRQ(ierr);
8407a35526eSPeter Brune     }
8417a35526eSPeter Brune   }
84294ae4db5SBarry Smith   ierr = PetscOptionsInt("-snes_linesearch_order","Order of approximation used in the line search","SNESLineSearchSetOrder",linesearch->order,&linesearch->order,NULL);CHKERRQ(ierr);
84394ae4db5SBarry Smith   ierr = PetscOptionsBool("-snes_linesearch_norms","Compute final norms in line search","SNESLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,NULL);CHKERRQ(ierr);
8447a35526eSPeter Brune 
8455a9b6599SPeter Brune   if (linesearch->ops->setfromoptions) {
846e55864a3SBarry Smith     (*linesearch->ops->setfromoptions)(PetscOptionsObject,linesearch);CHKERRQ(ierr);
8475a9b6599SPeter Brune   }
8485a9b6599SPeter Brune 
8490633abcbSJed Brown   ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)linesearch);CHKERRQ(ierr);
850bf7f4e0aSPeter Brune   ierr = PetscOptionsEnd();CHKERRQ(ierr);
851bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
852bf7f4e0aSPeter Brune }
853bf7f4e0aSPeter Brune 
854f40b411bSPeter Brune /*@
855f190f2fcSBarry Smith    SNESLineSearchView - Prints useful information about the line search
856f40b411bSPeter Brune 
857f40b411bSPeter Brune    Input Parameters:
8588e557f58SPeter Brune .  linesearch - linesearch context
859f40b411bSPeter Brune 
860f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
861f40b411bSPeter Brune 
862f40b411bSPeter Brune    Level: intermediate
863f40b411bSPeter Brune 
864dcf2fd19SBarry Smith .seealso: SNESLineSearchCreate()
865f40b411bSPeter Brune @*/
866bf388a1fSBarry Smith PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer)
867bf388a1fSBarry Smith {
8687f1410a3SPeter Brune   PetscErrorCode ierr;
8697f1410a3SPeter Brune   PetscBool      iascii;
870bf388a1fSBarry Smith 
871bf7f4e0aSPeter Brune   PetscFunctionBegin;
8727f1410a3SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
8737f1410a3SPeter Brune   if (!viewer) {
874ce94432eSBarry Smith     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch),&viewer);CHKERRQ(ierr);
8757f1410a3SPeter Brune   }
8767f1410a3SPeter Brune   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
8777f1410a3SPeter Brune   PetscCheckSameComm(linesearch,1,viewer,2);
878f40b411bSPeter Brune 
879251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
8807f1410a3SPeter Brune   if (iascii) {
881dae58748SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)linesearch,viewer);CHKERRQ(ierr);
8827f1410a3SPeter Brune     if (linesearch->ops->view) {
8837f1410a3SPeter Brune       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
8847f1410a3SPeter Brune       ierr = (*linesearch->ops->view)(linesearch,viewer);CHKERRQ(ierr);
8857f1410a3SPeter Brune       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
8867f1410a3SPeter Brune     }
887c69d1a72SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep,(double)linesearch->steptol);CHKERRQ(ierr);
888c69d1a72SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol,(double)linesearch->atol,(double)linesearch->ltol);CHKERRQ(ierr);
8897f1410a3SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D\n", linesearch->max_its);CHKERRQ(ierr);
8906b2b7091SBarry Smith     if (linesearch->ops->precheck) {
8916b2b7091SBarry Smith       if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) {
8927f1410a3SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  using precheck step to speed up Picard convergence\n", linesearch->max_its);CHKERRQ(ierr);
8937f1410a3SPeter Brune       } else {
8947f1410a3SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  using user-defined precheck step\n", linesearch->max_its);CHKERRQ(ierr);
8957f1410a3SPeter Brune       }
8967f1410a3SPeter Brune     }
8976b2b7091SBarry Smith     if (linesearch->ops->postcheck) {
8987f1410a3SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer,"  using user-defined postcheck step\n", linesearch->max_its);CHKERRQ(ierr);
8997f1410a3SPeter Brune     }
9007f1410a3SPeter Brune   }
901bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
902bf7f4e0aSPeter Brune }
903bf7f4e0aSPeter Brune 
904ea5d4fccSPeter Brune /*@C
905a80ff896SJed Brown    SNESLineSearchGetType - Gets the linesearch type
906a80ff896SJed Brown 
907a80ff896SJed Brown    Logically Collective on SNESLineSearch
908a80ff896SJed Brown 
909a80ff896SJed Brown    Input Parameters:
910a80ff896SJed Brown .  linesearch - linesearch context
911a80ff896SJed Brown 
912a80ff896SJed Brown    Output Parameters:
913a80ff896SJed Brown -  type - The type of line search, or NULL if not set
914a80ff896SJed Brown 
915a80ff896SJed Brown    Level: intermediate
916a80ff896SJed Brown 
917a80ff896SJed Brown .seealso: SNESLineSearchCreate(), SNESLineSearchType, SNESLineSearchSetFromOptions(), SNESLineSearchSetType()
918a80ff896SJed Brown @*/
919a80ff896SJed Brown PetscErrorCode SNESLineSearchGetType(SNESLineSearch linesearch, SNESLineSearchType *type)
920a80ff896SJed Brown {
921a80ff896SJed Brown   PetscFunctionBegin;
922a80ff896SJed Brown   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
923a80ff896SJed Brown   PetscValidCharPointer(type,2);
924a80ff896SJed Brown   *type = ((PetscObject)linesearch)->type_name;
925a80ff896SJed Brown   PetscFunctionReturn(0);
926a80ff896SJed Brown }
927a80ff896SJed Brown 
928a80ff896SJed Brown /*@C
929f1c6b773SPeter Brune    SNESLineSearchSetType - Sets the linesearch type
930f40b411bSPeter Brune 
931f190f2fcSBarry Smith    Logically Collective on SNESLineSearch
932f190f2fcSBarry Smith 
933f40b411bSPeter Brune    Input Parameters:
9348e557f58SPeter Brune +  linesearch - linesearch context
935f40b411bSPeter Brune -  type - The type of line search to be used
936f40b411bSPeter Brune 
937cd7522eaSPeter Brune    Available Types:
9381fe24845SBarry Smith +  SNESLINESEARCHBASIC - Simple damping line search, defaults to using the full Newton step
9391fe24845SBarry Smith .  SNESLINESEARCHBT - Backtracking line search over the L2 norm of the function
9401fe24845SBarry Smith .  SNESLINESEARCHL2 - Secant line search over the L2 norm of the function
9411fe24845SBarry Smith .  SNESLINESEARCHCP - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x)
9421fe24845SBarry Smith .  SNESLINESEARCHNLEQERR - Affine-covariant error-oriented linesearch
9431fe24845SBarry Smith -  SNESLINESEARCHSHELL - User provided SNESLineSearch implementation
9441fe24845SBarry Smith 
9451fe24845SBarry Smith    Options Database:
9461fe24845SBarry Smith .  -snes_linesearch_type <type> - basic, bt, l2, cp, nleqerr, shell
947cd7522eaSPeter Brune 
948f40b411bSPeter Brune    Level: intermediate
949f40b411bSPeter Brune 
950a80ff896SJed Brown .seealso: SNESLineSearchCreate(), SNESLineSearchType, SNESLineSearchSetFromOptions(), SNESLineSearchGetType()
951f40b411bSPeter Brune @*/
95219fd82e9SBarry Smith PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type)
953bf7f4e0aSPeter Brune {
954f1c6b773SPeter Brune   PetscErrorCode ierr,(*r)(SNESLineSearch);
955bf7f4e0aSPeter Brune   PetscBool      match;
956bf7f4e0aSPeter Brune 
957bf7f4e0aSPeter Brune   PetscFunctionBegin;
958f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
959bf7f4e0aSPeter Brune   PetscValidCharPointer(type,2);
960bf7f4e0aSPeter Brune 
961251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)linesearch,type,&match);CHKERRQ(ierr);
962bf7f4e0aSPeter Brune   if (match) PetscFunctionReturn(0);
963bf7f4e0aSPeter Brune 
9641c9cd337SJed Brown   ierr = PetscFunctionListFind(SNESLineSearchList,type,&r);CHKERRQ(ierr);
965bf7f4e0aSPeter Brune   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type);
966bf7f4e0aSPeter Brune   /* Destroy the previous private linesearch context */
967bf7f4e0aSPeter Brune   if (linesearch->ops->destroy) {
968bf7f4e0aSPeter Brune     ierr = (*(linesearch)->ops->destroy)(linesearch);CHKERRQ(ierr);
9690298fd71SBarry Smith     linesearch->ops->destroy = NULL;
970bf7f4e0aSPeter Brune   }
971f1c6b773SPeter Brune   /* Reinitialize function pointers in SNESLineSearchOps structure */
9729e5d0892SLisandro Dalcin   linesearch->ops->apply          = NULL;
9739e5d0892SLisandro Dalcin   linesearch->ops->view           = NULL;
9749e5d0892SLisandro Dalcin   linesearch->ops->setfromoptions = NULL;
9759e5d0892SLisandro Dalcin   linesearch->ops->destroy        = NULL;
976bf7f4e0aSPeter Brune 
977bf7f4e0aSPeter Brune   ierr = PetscObjectChangeTypeName((PetscObject)linesearch,type);CHKERRQ(ierr);
978bf7f4e0aSPeter Brune   ierr = (*r)(linesearch);CHKERRQ(ierr);
979bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
980bf7f4e0aSPeter Brune }
981bf7f4e0aSPeter Brune 
982f40b411bSPeter Brune /*@
98378bcb3b5SPeter Brune    SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation.
984f40b411bSPeter Brune 
985f40b411bSPeter Brune    Input Parameters:
9868e557f58SPeter Brune +  linesearch - linesearch context
987f40b411bSPeter Brune -  snes - The snes instance
988f40b411bSPeter Brune 
98978bcb3b5SPeter Brune    Level: developer
99078bcb3b5SPeter Brune 
99178bcb3b5SPeter Brune    Notes:
992f190f2fcSBarry Smith    This happens automatically when the line search is obtained/created with
9937601faf0SJed Brown    SNESGetLineSearch().  This routine is therefore mainly called within SNES
99478bcb3b5SPeter Brune    implementations.
995f40b411bSPeter Brune 
9968141a3b9SPeter Brune    Level: developer
997f40b411bSPeter Brune 
998cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
999f40b411bSPeter Brune @*/
1000bf388a1fSBarry Smith PetscErrorCode  SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes)
1001bf388a1fSBarry Smith {
1002bf7f4e0aSPeter Brune   PetscFunctionBegin;
1003f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1004bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,2);
1005bf7f4e0aSPeter Brune   linesearch->snes = snes;
1006bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1007bf7f4e0aSPeter Brune }
1008bf7f4e0aSPeter Brune 
1009f40b411bSPeter Brune /*@
10108141a3b9SPeter Brune    SNESLineSearchGetSNES - Gets the SNES instance associated with the line search.
10118141a3b9SPeter Brune    Having an associated SNES is necessary because most line search implementations must be able to
10128141a3b9SPeter Brune    evaluate the function using SNESComputeFunction() for the associated SNES.  This routine
10138141a3b9SPeter Brune    is used in the line search implementations when one must get this associated SNES instance.
1014f40b411bSPeter Brune 
1015f40b411bSPeter Brune    Input Parameters:
10168e557f58SPeter Brune .  linesearch - linesearch context
1017f40b411bSPeter Brune 
1018f40b411bSPeter Brune    Output Parameters:
1019f40b411bSPeter Brune .  snes - The snes instance
1020f40b411bSPeter Brune 
10218141a3b9SPeter Brune    Level: developer
1022f40b411bSPeter Brune 
1023cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
1024f40b411bSPeter Brune @*/
1025bf388a1fSBarry Smith PetscErrorCode  SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes)
1026bf388a1fSBarry Smith {
1027bf7f4e0aSPeter Brune   PetscFunctionBegin;
1028f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
10296a388c36SPeter Brune   PetscValidPointer(snes, 2);
1030bf7f4e0aSPeter Brune   *snes = linesearch->snes;
1031bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1032bf7f4e0aSPeter Brune }
1033bf7f4e0aSPeter Brune 
1034f40b411bSPeter Brune /*@
1035f1c6b773SPeter Brune    SNESLineSearchGetLambda - Gets the last linesearch steplength discovered.
1036f40b411bSPeter Brune 
1037f40b411bSPeter Brune    Input Parameters:
10388e557f58SPeter Brune .  linesearch - linesearch context
1039f40b411bSPeter Brune 
1040f40b411bSPeter Brune    Output Parameters:
1041cd7522eaSPeter Brune .  lambda - The last steplength computed during SNESLineSearchApply()
1042f40b411bSPeter Brune 
104378bcb3b5SPeter Brune    Level: advanced
104478bcb3b5SPeter Brune 
10458e557f58SPeter Brune    Notes:
10468e557f58SPeter Brune    This is useful in methods where the solver is ill-scaled and
104778bcb3b5SPeter Brune    requires some adaptive notion of the difference in scale between the
104878bcb3b5SPeter Brune    solution and the function.  For instance, SNESQN may be scaled by the
104978bcb3b5SPeter Brune    line search lambda using the argument -snes_qn_scaling ls.
105078bcb3b5SPeter Brune 
1051cd7522eaSPeter Brune .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetDamping(), SNESLineSearchApply()
1052f40b411bSPeter Brune @*/
1053f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda)
10546a388c36SPeter Brune {
10556a388c36SPeter Brune   PetscFunctionBegin;
1056f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1057534a8f05SLisandro Dalcin   PetscValidRealPointer(lambda, 2);
10586a388c36SPeter Brune   *lambda = linesearch->lambda;
10596a388c36SPeter Brune   PetscFunctionReturn(0);
10606a388c36SPeter Brune }
10616a388c36SPeter Brune 
1062f40b411bSPeter Brune /*@
1063f1c6b773SPeter Brune    SNESLineSearchSetLambda - Sets the linesearch steplength.
1064f40b411bSPeter Brune 
1065f40b411bSPeter Brune    Input Parameters:
10668e557f58SPeter Brune +  linesearch - linesearch context
1067f40b411bSPeter Brune -  lambda - The last steplength.
1068f40b411bSPeter Brune 
1069cd7522eaSPeter Brune    Notes:
1070f190f2fcSBarry Smith    This routine is typically used within implementations of SNESLineSearchApply()
1071cd7522eaSPeter Brune    to set the final steplength.  This routine (and SNESLineSearchGetLambda()) were
1072cd7522eaSPeter Brune    added in order to facilitate Quasi-Newton methods that use the previous steplength
1073cd7522eaSPeter Brune    as an inner scaling parameter.
1074cd7522eaSPeter Brune 
107578bcb3b5SPeter Brune    Level: advanced
1076f40b411bSPeter Brune 
1077f1c6b773SPeter Brune .seealso: SNESLineSearchGetLambda()
1078f40b411bSPeter Brune @*/
1079f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
10806a388c36SPeter Brune {
10816a388c36SPeter Brune   PetscFunctionBegin;
1082f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
10836a388c36SPeter Brune   linesearch->lambda = lambda;
10846a388c36SPeter Brune   PetscFunctionReturn(0);
10856a388c36SPeter Brune }
10866a388c36SPeter Brune 
1087f40b411bSPeter Brune /*@
10883c7d6663SPeter Brune    SNESLineSearchGetTolerances - Gets the tolerances for the linesearch.  These include
108978bcb3b5SPeter Brune    tolerances for the relative and absolute change in the function norm, the change
109078bcb3b5SPeter Brune    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
109178bcb3b5SPeter Brune    and the maximum number of iterations the line search procedure may take.
1092f40b411bSPeter Brune 
1093f40b411bSPeter Brune    Input Parameters:
10948e557f58SPeter Brune .  linesearch - linesearch context
1095f40b411bSPeter Brune 
1096f40b411bSPeter Brune    Output Parameters:
1097516fe3c3SPeter Brune +  steptol - The minimum steplength
10986cc8e53bSPeter Brune .  maxstep - The maximum steplength
1099516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
1100516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
1101516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
1102516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
1103f40b411bSPeter Brune 
110478bcb3b5SPeter Brune    Level: intermediate
110578bcb3b5SPeter Brune 
110678bcb3b5SPeter Brune    Notes:
110778bcb3b5SPeter Brune    Different line searches may implement these parameters slightly differently as
11083c7d6663SPeter Brune    the type requires.
1109516fe3c3SPeter Brune 
1110f1c6b773SPeter Brune .seealso: SNESLineSearchSetTolerances()
1111f40b411bSPeter Brune @*/
1112f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
11136a388c36SPeter Brune {
11146a388c36SPeter Brune   PetscFunctionBegin;
1115f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1116516fe3c3SPeter Brune   if (steptol) {
1117534a8f05SLisandro Dalcin     PetscValidRealPointer(steptol, 2);
11186a388c36SPeter Brune     *steptol = linesearch->steptol;
1119516fe3c3SPeter Brune   }
1120516fe3c3SPeter Brune   if (maxstep) {
1121534a8f05SLisandro Dalcin     PetscValidRealPointer(maxstep, 3);
1122516fe3c3SPeter Brune     *maxstep = linesearch->maxstep;
1123516fe3c3SPeter Brune   }
1124516fe3c3SPeter Brune   if (rtol) {
1125534a8f05SLisandro Dalcin     PetscValidRealPointer(rtol, 4);
1126516fe3c3SPeter Brune     *rtol = linesearch->rtol;
1127516fe3c3SPeter Brune   }
1128516fe3c3SPeter Brune   if (atol) {
1129534a8f05SLisandro Dalcin     PetscValidRealPointer(atol, 5);
1130516fe3c3SPeter Brune     *atol = linesearch->atol;
1131516fe3c3SPeter Brune   }
1132516fe3c3SPeter Brune   if (ltol) {
1133534a8f05SLisandro Dalcin     PetscValidRealPointer(ltol, 6);
1134516fe3c3SPeter Brune     *ltol = linesearch->ltol;
1135516fe3c3SPeter Brune   }
1136516fe3c3SPeter Brune   if (max_its) {
1137534a8f05SLisandro Dalcin     PetscValidIntPointer(max_its, 7);
1138516fe3c3SPeter Brune     *max_its = linesearch->max_its;
1139516fe3c3SPeter Brune   }
11406a388c36SPeter Brune   PetscFunctionReturn(0);
11416a388c36SPeter Brune }
11426a388c36SPeter Brune 
1143f40b411bSPeter Brune /*@
11443c7d6663SPeter Brune    SNESLineSearchSetTolerances -  Gets the tolerances for the linesearch.  These include
114578bcb3b5SPeter Brune    tolerances for the relative and absolute change in the function norm, the change
114678bcb3b5SPeter Brune    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
114778bcb3b5SPeter Brune    and the maximum number of iterations the line search procedure may take.
1148f40b411bSPeter Brune 
1149f40b411bSPeter Brune    Input Parameters:
11508e557f58SPeter Brune +  linesearch - linesearch context
1151516fe3c3SPeter Brune .  steptol - The minimum steplength
11526cc8e53bSPeter Brune .  maxstep - The maximum steplength
1153516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
1154516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
1155516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
1156516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
1157f40b411bSPeter Brune 
115878bcb3b5SPeter Brune    Notes:
11593c7d6663SPeter Brune    The user may choose to not set any of the tolerances using PETSC_DEFAULT in
116078bcb3b5SPeter Brune    place of an argument.
1161f40b411bSPeter Brune 
116278bcb3b5SPeter Brune    Level: intermediate
1163516fe3c3SPeter Brune 
1164f1c6b773SPeter Brune .seealso: SNESLineSearchGetTolerances()
1165f40b411bSPeter Brune @*/
1166f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its)
11676a388c36SPeter Brune {
11686a388c36SPeter Brune   PetscFunctionBegin;
1169f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1170d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,steptol,2);
1171d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,maxstep,3);
1172d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,rtol,4);
1173d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,atol,5);
1174d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,ltol,6);
1175d3952184SSatish Balay   PetscValidLogicalCollectiveInt(linesearch,max_its,7);
1176d3952184SSatish Balay 
1177d3952184SSatish Balay   if (steptol!= PETSC_DEFAULT) {
1178ce94432eSBarry Smith     if (steptol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Minimum step length %14.12e must be non-negative",(double)steptol);
11796a388c36SPeter Brune     linesearch->steptol = steptol;
1180d3952184SSatish Balay   }
1181d3952184SSatish Balay 
1182d3952184SSatish Balay   if (maxstep!= PETSC_DEFAULT) {
1183ce94432eSBarry Smith     if (maxstep < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum step length %14.12e must be non-negative",(double)maxstep);
1184516fe3c3SPeter Brune     linesearch->maxstep = maxstep;
1185d3952184SSatish Balay   }
1186d3952184SSatish Balay 
1187d3952184SSatish Balay   if (rtol != PETSC_DEFAULT) {
1188ce94432eSBarry 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);
1189516fe3c3SPeter Brune     linesearch->rtol = rtol;
1190d3952184SSatish Balay   }
1191d3952184SSatish Balay 
1192d3952184SSatish Balay   if (atol != PETSC_DEFAULT) {
1193ce94432eSBarry Smith     if (atol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %14.12e must be non-negative",(double)atol);
1194516fe3c3SPeter Brune     linesearch->atol = atol;
1195d3952184SSatish Balay   }
1196d3952184SSatish Balay 
1197d3952184SSatish Balay   if (ltol != PETSC_DEFAULT) {
11982479783cSJose E. Roman     if (ltol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Lambda tolerance %14.12e must be non-negative",(double)ltol);
1199516fe3c3SPeter Brune     linesearch->ltol = ltol;
1200d3952184SSatish Balay   }
1201d3952184SSatish Balay 
1202d3952184SSatish Balay   if (max_its != PETSC_DEFAULT) {
1203ce94432eSBarry Smith     if (max_its < 0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",max_its);
1204516fe3c3SPeter Brune     linesearch->max_its = max_its;
1205d3952184SSatish Balay   }
12066a388c36SPeter Brune   PetscFunctionReturn(0);
12076a388c36SPeter Brune }
12086a388c36SPeter Brune 
1209f40b411bSPeter Brune /*@
1210f1c6b773SPeter Brune    SNESLineSearchGetDamping - Gets the line search damping parameter.
1211f40b411bSPeter Brune 
1212f40b411bSPeter Brune    Input Parameters:
12138e557f58SPeter Brune .  linesearch - linesearch context
1214f40b411bSPeter Brune 
1215f40b411bSPeter Brune    Output Parameters:
12168e557f58SPeter Brune .  damping - The damping parameter
1217f40b411bSPeter Brune 
121878bcb3b5SPeter Brune    Level: advanced
1219f40b411bSPeter Brune 
122078bcb3b5SPeter Brune .seealso: SNESLineSearchGetStepTolerance(), SNESQN
1221f40b411bSPeter Brune @*/
1222f40b411bSPeter Brune 
1223f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping)
12246a388c36SPeter Brune {
12256a388c36SPeter Brune   PetscFunctionBegin;
1226f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1227534a8f05SLisandro Dalcin   PetscValidRealPointer(damping, 2);
12286a388c36SPeter Brune   *damping = linesearch->damping;
12296a388c36SPeter Brune   PetscFunctionReturn(0);
12306a388c36SPeter Brune }
12316a388c36SPeter Brune 
1232f40b411bSPeter Brune /*@
1233fd292e60Sprj-    SNESLineSearchSetDamping - Sets the line search damping parameter.
1234f40b411bSPeter Brune 
1235f40b411bSPeter Brune    Input Parameters:
123603fd4120SBarry Smith +  linesearch - linesearch context
123703fd4120SBarry Smith -  damping - The damping parameter
1238f40b411bSPeter Brune 
123903fd4120SBarry Smith    Options Database:
124003fd4120SBarry Smith .   -snes_linesearch_damping
1241f40b411bSPeter Brune    Level: intermediate
1242f40b411bSPeter Brune 
1243cd7522eaSPeter Brune    Notes:
1244cd7522eaSPeter Brune    The basic line search merely takes the update step scaled by the damping parameter.
1245cd7522eaSPeter Brune    The use of the damping parameter in the l2 and cp line searches is much more subtle;
124678bcb3b5SPeter Brune    it is used as a starting point in calculating the secant step. However, the eventual
1247cd7522eaSPeter Brune    step may be of greater length than the damping parameter.  In the bt line search it is
1248cd7522eaSPeter Brune    used as the maximum possible step length, as the bt line search only backtracks.
1249cd7522eaSPeter Brune 
1250f1c6b773SPeter Brune .seealso: SNESLineSearchGetDamping()
1251f40b411bSPeter Brune @*/
1252f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping)
12536a388c36SPeter Brune {
12546a388c36SPeter Brune   PetscFunctionBegin;
1255f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
12566a388c36SPeter Brune   linesearch->damping = damping;
12576a388c36SPeter Brune   PetscFunctionReturn(0);
12586a388c36SPeter Brune }
12596a388c36SPeter Brune 
126059405d5eSPeter Brune /*@
126159405d5eSPeter Brune    SNESLineSearchGetOrder - Gets the line search approximation order.
126259405d5eSPeter Brune 
126359405d5eSPeter Brune    Input Parameters:
126478bcb3b5SPeter Brune .  linesearch - linesearch context
126559405d5eSPeter Brune 
126659405d5eSPeter Brune    Output Parameters:
12678e557f58SPeter Brune .  order - The order
126859405d5eSPeter Brune 
126978bcb3b5SPeter Brune    Possible Values for order:
12703c7d6663SPeter Brune +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
12713c7d6663SPeter Brune .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
12723c7d6663SPeter Brune -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
127378bcb3b5SPeter Brune 
127459405d5eSPeter Brune    Level: intermediate
127559405d5eSPeter Brune 
127659405d5eSPeter Brune .seealso: SNESLineSearchSetOrder()
127759405d5eSPeter Brune @*/
127859405d5eSPeter Brune 
1279b000cd8dSPeter Brune PetscErrorCode  SNESLineSearchGetOrder(SNESLineSearch linesearch,PetscInt *order)
128059405d5eSPeter Brune {
128159405d5eSPeter Brune   PetscFunctionBegin;
128259405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1283534a8f05SLisandro Dalcin   PetscValidIntPointer(order, 2);
128459405d5eSPeter Brune   *order = linesearch->order;
128559405d5eSPeter Brune   PetscFunctionReturn(0);
128659405d5eSPeter Brune }
128759405d5eSPeter Brune 
128859405d5eSPeter Brune /*@
12891f8196a2SJed Brown    SNESLineSearchSetOrder - Sets the maximum order of the polynomial fit used in the line search
129059405d5eSPeter Brune 
129159405d5eSPeter Brune    Input Parameters:
1292a2b725a8SWilliam Gropp +  linesearch - linesearch context
1293a2b725a8SWilliam Gropp -  order - The damping parameter
129459405d5eSPeter Brune 
129559405d5eSPeter Brune    Level: intermediate
129659405d5eSPeter Brune 
129778bcb3b5SPeter Brune    Possible Values for order:
12983c7d6663SPeter Brune +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
12993c7d6663SPeter Brune .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
13003c7d6663SPeter Brune -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
130178bcb3b5SPeter Brune 
130259405d5eSPeter Brune    Notes:
130359405d5eSPeter Brune    Variable orders are supported by the following line searches:
130478bcb3b5SPeter Brune +  bt - cubic and quadratic
130578bcb3b5SPeter Brune -  cp - linear and quadratic
130659405d5eSPeter Brune 
13071f8196a2SJed Brown .seealso: SNESLineSearchGetOrder(), SNESLineSearchSetDamping()
130859405d5eSPeter Brune @*/
1309b000cd8dSPeter Brune PetscErrorCode  SNESLineSearchSetOrder(SNESLineSearch linesearch,PetscInt order)
131059405d5eSPeter Brune {
131159405d5eSPeter Brune   PetscFunctionBegin;
131259405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
131359405d5eSPeter Brune   linesearch->order = order;
131459405d5eSPeter Brune   PetscFunctionReturn(0);
131559405d5eSPeter Brune }
131659405d5eSPeter Brune 
1317f40b411bSPeter Brune /*@
1318f1c6b773SPeter Brune    SNESLineSearchGetNorms - Gets the norms for for X, Y, and F.
1319f40b411bSPeter Brune 
1320f40b411bSPeter Brune    Input Parameters:
132178bcb3b5SPeter Brune .  linesearch - linesearch context
1322f40b411bSPeter Brune 
1323f40b411bSPeter Brune    Output Parameters:
1324f40b411bSPeter Brune +  xnorm - The norm of the current solution
1325f40b411bSPeter Brune .  fnorm - The norm of the current function
1326f40b411bSPeter Brune -  ynorm - The norm of the current update
1327f40b411bSPeter Brune 
1328cd7522eaSPeter Brune    Notes:
1329cd7522eaSPeter Brune    This function is mainly called from SNES implementations.
1330cd7522eaSPeter Brune 
133178bcb3b5SPeter Brune    Level: developer
1332f40b411bSPeter Brune 
1333f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms() SNESLineSearchGetVecs()
1334f40b411bSPeter Brune @*/
1335f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm)
1336bf7f4e0aSPeter Brune {
1337bf7f4e0aSPeter Brune   PetscFunctionBegin;
1338f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1339f5af7f23SKarl Rupp   if (xnorm) *xnorm = linesearch->xnorm;
1340f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
1341f5af7f23SKarl Rupp   if (ynorm) *ynorm = linesearch->ynorm;
1342bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1343bf7f4e0aSPeter Brune }
1344bf7f4e0aSPeter Brune 
1345f40b411bSPeter Brune /*@
1346f1c6b773SPeter Brune    SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F.
1347f40b411bSPeter Brune 
1348f40b411bSPeter Brune    Input Parameters:
134978bcb3b5SPeter Brune +  linesearch - linesearch context
1350f40b411bSPeter Brune .  xnorm - The norm of the current solution
1351f40b411bSPeter Brune .  fnorm - The norm of the current function
1352f40b411bSPeter Brune -  ynorm - The norm of the current update
1353f40b411bSPeter Brune 
135478bcb3b5SPeter Brune    Level: advanced
1355f40b411bSPeter Brune 
1356f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1357f40b411bSPeter Brune @*/
1358f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
13596a388c36SPeter Brune {
13606a388c36SPeter Brune   PetscFunctionBegin;
1361f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
13626a388c36SPeter Brune   linesearch->xnorm = xnorm;
13636a388c36SPeter Brune   linesearch->fnorm = fnorm;
13646a388c36SPeter Brune   linesearch->ynorm = ynorm;
13656a388c36SPeter Brune   PetscFunctionReturn(0);
13666a388c36SPeter Brune }
13676a388c36SPeter Brune 
1368f40b411bSPeter Brune /*@
1369f1c6b773SPeter Brune    SNESLineSearchComputeNorms - Computes the norms of X, F, and Y.
1370f40b411bSPeter Brune 
1371f40b411bSPeter Brune    Input Parameters:
137278bcb3b5SPeter Brune .  linesearch - linesearch context
1373f40b411bSPeter Brune 
1374f40b411bSPeter Brune    Options Database Keys:
13758e557f58SPeter Brune .   -snes_linesearch_norms - turn norm computation on or off
1376f40b411bSPeter Brune 
1377f40b411bSPeter Brune    Level: intermediate
1378f40b411bSPeter Brune 
137978bcb3b5SPeter Brune .seealso: SNESLineSearchGetNorms, SNESLineSearchSetNorms(), SNESLineSearchSetComputeNorms()
1380f40b411bSPeter Brune @*/
1381f1c6b773SPeter Brune PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
13826a388c36SPeter Brune {
13836a388c36SPeter Brune   PetscErrorCode ierr;
13849bd66eb0SPeter Brune   SNES           snes;
1385bf388a1fSBarry Smith 
13866a388c36SPeter Brune   PetscFunctionBegin;
13876a388c36SPeter Brune   if (linesearch->norms) {
13889bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
1389f1c6b773SPeter Brune       ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
13909bd66eb0SPeter Brune       ierr = VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
13919bd66eb0SPeter Brune       ierr = VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
13929bd66eb0SPeter Brune       ierr = (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);CHKERRQ(ierr);
13939bd66eb0SPeter Brune     } else {
13946a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_func,   NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
13956a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_sol,    NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
13966a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
13976a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_func,     NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
13986a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_sol,      NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
13996a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_update,   NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
14006a388c36SPeter Brune     }
14019bd66eb0SPeter Brune   }
14026a388c36SPeter Brune   PetscFunctionReturn(0);
14036a388c36SPeter Brune }
14046a388c36SPeter Brune 
14056f263ca3SPeter Brune /*@
14066f263ca3SPeter Brune    SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search.
14076f263ca3SPeter Brune 
14086f263ca3SPeter Brune    Input Parameters:
140978bcb3b5SPeter Brune +  linesearch  - linesearch context
141078bcb3b5SPeter Brune -  flg  - indicates whether or not to compute norms
14116f263ca3SPeter Brune 
14126f263ca3SPeter Brune    Options Database Keys:
14131f8196a2SJed Brown .   -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic linesearch
14146f263ca3SPeter Brune 
14156f263ca3SPeter Brune    Notes:
14161f8196a2SJed Brown    This is most relevant to the SNESLINESEARCHBASIC line search type since most line searches have a stopping criteria involving the norm.
14176f263ca3SPeter Brune 
14186f263ca3SPeter Brune    Level: intermediate
14196f263ca3SPeter Brune 
14201a4f838cSPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetNorms(), SNESLineSearchComputeNorms(), SNESLINESEARCHBASIC
14216f263ca3SPeter Brune @*/
14226f263ca3SPeter Brune PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg)
14236f263ca3SPeter Brune {
14246f263ca3SPeter Brune   PetscFunctionBegin;
14256f263ca3SPeter Brune   linesearch->norms = flg;
14266f263ca3SPeter Brune   PetscFunctionReturn(0);
14276f263ca3SPeter Brune }
14286f263ca3SPeter Brune 
1429f40b411bSPeter Brune /*@
1430f1c6b773SPeter Brune    SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context
1431f40b411bSPeter Brune 
1432f40b411bSPeter Brune    Input Parameters:
143378bcb3b5SPeter Brune .  linesearch - linesearch context
1434f40b411bSPeter Brune 
1435f40b411bSPeter Brune    Output Parameters:
14366232e825SPeter Brune +  X - Solution vector
14376232e825SPeter Brune .  F - Function vector
14386232e825SPeter Brune .  Y - Search direction vector
14396232e825SPeter Brune .  W - Solution work vector
14406232e825SPeter Brune -  G - Function work vector
14416232e825SPeter Brune 
14427bba9028SPeter Brune    Notes:
14437bba9028SPeter Brune    At the beginning of a line search application, X should contain a
14446232e825SPeter Brune    solution and the vector F the function computed at X.  At the end of the
14456232e825SPeter Brune    line search application, X should contain the new solution, and F the
14466232e825SPeter Brune    function evaluated at the new solution.
1447f40b411bSPeter Brune 
14482a7a6963SBarry Smith    These vectors are owned by the SNESLineSearch and should not be destroyed by the caller
14492a7a6963SBarry Smith 
145078bcb3b5SPeter Brune    Level: advanced
1451f40b411bSPeter Brune 
1452f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1453f40b411bSPeter Brune @*/
1454bf388a1fSBarry Smith PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G)
1455bf388a1fSBarry Smith {
14566a388c36SPeter Brune   PetscFunctionBegin;
1457f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
14586a388c36SPeter Brune   if (X) {
14596a388c36SPeter Brune     PetscValidPointer(X, 2);
14606a388c36SPeter Brune     *X = linesearch->vec_sol;
14616a388c36SPeter Brune   }
14626a388c36SPeter Brune   if (F) {
14636a388c36SPeter Brune     PetscValidPointer(F, 3);
14646a388c36SPeter Brune     *F = linesearch->vec_func;
14656a388c36SPeter Brune   }
14666a388c36SPeter Brune   if (Y) {
14676a388c36SPeter Brune     PetscValidPointer(Y, 4);
14686a388c36SPeter Brune     *Y = linesearch->vec_update;
14696a388c36SPeter Brune   }
14706a388c36SPeter Brune   if (W) {
14716a388c36SPeter Brune     PetscValidPointer(W, 5);
14726a388c36SPeter Brune     *W = linesearch->vec_sol_new;
14736a388c36SPeter Brune   }
14746a388c36SPeter Brune   if (G) {
14756a388c36SPeter Brune     PetscValidPointer(G, 6);
14766a388c36SPeter Brune     *G = linesearch->vec_func_new;
14776a388c36SPeter Brune   }
14786a388c36SPeter Brune   PetscFunctionReturn(0);
14796a388c36SPeter Brune }
14806a388c36SPeter Brune 
1481f40b411bSPeter Brune /*@
1482f1c6b773SPeter Brune    SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context
1483f40b411bSPeter Brune 
1484f40b411bSPeter Brune    Input Parameters:
148578bcb3b5SPeter Brune +  linesearch - linesearch context
14866232e825SPeter Brune .  X - Solution vector
14876232e825SPeter Brune .  F - Function vector
14886232e825SPeter Brune .  Y - Search direction vector
14896232e825SPeter Brune .  W - Solution work vector
14906232e825SPeter Brune -  G - Function work vector
1491f40b411bSPeter Brune 
149278bcb3b5SPeter Brune    Level: advanced
1493f40b411bSPeter Brune 
1494f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms(), SNESLineSearchGetVecs()
1495f40b411bSPeter Brune @*/
1496bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G)
1497bf388a1fSBarry Smith {
14986a388c36SPeter Brune   PetscFunctionBegin;
1499f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
15006a388c36SPeter Brune   if (X) {
15016a388c36SPeter Brune     PetscValidHeaderSpecific(X,VEC_CLASSID,2);
15026a388c36SPeter Brune     linesearch->vec_sol = X;
15036a388c36SPeter Brune   }
15046a388c36SPeter Brune   if (F) {
15056a388c36SPeter Brune     PetscValidHeaderSpecific(F,VEC_CLASSID,3);
15066a388c36SPeter Brune     linesearch->vec_func = F;
15076a388c36SPeter Brune   }
15086a388c36SPeter Brune   if (Y) {
15096a388c36SPeter Brune     PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
15106a388c36SPeter Brune     linesearch->vec_update = Y;
15116a388c36SPeter Brune   }
15126a388c36SPeter Brune   if (W) {
15136a388c36SPeter Brune     PetscValidHeaderSpecific(W,VEC_CLASSID,5);
15146a388c36SPeter Brune     linesearch->vec_sol_new = W;
15156a388c36SPeter Brune   }
15166a388c36SPeter Brune   if (G) {
15176a388c36SPeter Brune     PetscValidHeaderSpecific(G,VEC_CLASSID,6);
15186a388c36SPeter Brune     linesearch->vec_func_new = G;
15196a388c36SPeter Brune   }
15206a388c36SPeter Brune   PetscFunctionReturn(0);
15216a388c36SPeter Brune }
15226a388c36SPeter Brune 
1523e7058c64SPeter Brune /*@C
1524f1c6b773SPeter Brune    SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1525e7058c64SPeter Brune    SNES options in the database.
1526e7058c64SPeter Brune 
1527cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
1528e7058c64SPeter Brune 
1529e7058c64SPeter Brune    Input Parameters:
1530e7058c64SPeter Brune +  snes - the SNES context
1531e7058c64SPeter Brune -  prefix - the prefix to prepend to all option names
1532e7058c64SPeter Brune 
1533e7058c64SPeter Brune    Notes:
1534e7058c64SPeter Brune    A hyphen (-) must NOT be given at the beginning of the prefix name.
1535e7058c64SPeter Brune    The first character of all runtime options is AUTOMATICALLY the hyphen.
1536e7058c64SPeter Brune 
1537e7058c64SPeter Brune    Level: advanced
1538e7058c64SPeter Brune 
1539e7058c64SPeter Brune .seealso: SNESGetOptionsPrefix()
1540e7058c64SPeter Brune @*/
1541f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[])
1542e7058c64SPeter Brune {
1543e7058c64SPeter Brune   PetscErrorCode ierr;
1544e7058c64SPeter Brune 
1545e7058c64SPeter Brune   PetscFunctionBegin;
1546f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1547e7058c64SPeter Brune   ierr = PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1548e7058c64SPeter Brune   PetscFunctionReturn(0);
1549e7058c64SPeter Brune }
1550e7058c64SPeter Brune 
1551e7058c64SPeter Brune /*@C
1552f1c6b773SPeter Brune    SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all
1553f1c6b773SPeter Brune    SNESLineSearch options in the database.
1554e7058c64SPeter Brune 
1555e7058c64SPeter Brune    Not Collective
1556e7058c64SPeter Brune 
1557e7058c64SPeter Brune    Input Parameter:
1558cd7522eaSPeter Brune .  linesearch - the SNESLineSearch context
1559e7058c64SPeter Brune 
1560e7058c64SPeter Brune    Output Parameter:
1561e7058c64SPeter Brune .  prefix - pointer to the prefix string used
1562e7058c64SPeter Brune 
15638e557f58SPeter Brune    Notes:
15648e557f58SPeter Brune    On the fortran side, the user should pass in a string 'prefix' of
1565e7058c64SPeter Brune    sufficient length to hold the prefix.
1566e7058c64SPeter Brune 
1567e7058c64SPeter Brune    Level: advanced
1568e7058c64SPeter Brune 
1569e7058c64SPeter Brune .seealso: SNESAppendOptionsPrefix()
1570e7058c64SPeter Brune @*/
1571f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[])
1572e7058c64SPeter Brune {
1573e7058c64SPeter Brune   PetscErrorCode ierr;
1574e7058c64SPeter Brune 
1575e7058c64SPeter Brune   PetscFunctionBegin;
1576f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1577e7058c64SPeter Brune   ierr = PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1578e7058c64SPeter Brune   PetscFunctionReturn(0);
1579e7058c64SPeter Brune }
1580bf7f4e0aSPeter Brune 
15818d359177SBarry Smith /*@C
1582fa0ddf94SBarry Smith    SNESLineSearchSetWorkVecs - Gets work vectors for the line search.
1583f40b411bSPeter Brune 
1584*d8d19677SJose E. Roman    Input Parameters:
1585f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
1586f40b411bSPeter Brune -  nwork - the number of work vectors
1587f40b411bSPeter Brune 
1588f40b411bSPeter Brune    Level: developer
1589f40b411bSPeter Brune 
1590fa0ddf94SBarry Smith .seealso: SNESSetWorkVecs()
1591f40b411bSPeter Brune @*/
1592fa0ddf94SBarry Smith PetscErrorCode  SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork)
1593bf7f4e0aSPeter Brune {
1594bf7f4e0aSPeter Brune   PetscErrorCode ierr;
1595bf388a1fSBarry Smith 
1596bf7f4e0aSPeter Brune   PetscFunctionBegin;
1597bf7f4e0aSPeter Brune   if (linesearch->vec_sol) {
1598bf7f4e0aSPeter Brune     ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr);
15998d359177SBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
1600bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1601bf7f4e0aSPeter Brune }
1602bf7f4e0aSPeter Brune 
1603f40b411bSPeter Brune /*@
1604422a814eSBarry Smith    SNESLineSearchGetReason - Gets the success/failure status of the last line search application
1605f40b411bSPeter Brune 
1606f40b411bSPeter Brune    Input Parameters:
160778bcb3b5SPeter Brune .  linesearch - linesearch context
1608f40b411bSPeter Brune 
1609f40b411bSPeter Brune    Output Parameters:
1610422a814eSBarry Smith .  result - The success or failure status
1611f40b411bSPeter Brune 
1612cd7522eaSPeter Brune    Notes:
1613c98378a5SBarry Smith    This is typically called after SNESLineSearchApply() in order to determine if the line-search failed
1614cd7522eaSPeter Brune    (and set the SNES convergence accordingly).
1615cd7522eaSPeter Brune 
1616f40b411bSPeter Brune    Level: intermediate
1617f40b411bSPeter Brune 
1618422a814eSBarry Smith .seealso: SNESLineSearchSetReason(), SNESLineSearchReason
1619f40b411bSPeter Brune @*/
1620422a814eSBarry Smith PetscErrorCode  SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result)
1621bf7f4e0aSPeter Brune {
1622bf7f4e0aSPeter Brune   PetscFunctionBegin;
1623f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1624422a814eSBarry Smith   PetscValidPointer(result, 2);
1625422a814eSBarry Smith   *result = linesearch->result;
1626bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1627bf7f4e0aSPeter Brune }
1628bf7f4e0aSPeter Brune 
1629f40b411bSPeter Brune /*@
1630422a814eSBarry Smith    SNESLineSearchSetReason - Sets the success/failure status of the last line search application
1631f40b411bSPeter Brune 
1632f40b411bSPeter Brune    Input Parameters:
163378bcb3b5SPeter Brune +  linesearch - linesearch context
1634422a814eSBarry Smith -  result - The success or failure status
1635f40b411bSPeter Brune 
1636cd7522eaSPeter Brune    Notes:
1637cd7522eaSPeter Brune    This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set
1638cd7522eaSPeter Brune    the success or failure of the line search method.
1639cd7522eaSPeter Brune 
164078bcb3b5SPeter Brune    Level: developer
1641f40b411bSPeter Brune 
1642422a814eSBarry Smith .seealso: SNESLineSearchGetSResult()
1643f40b411bSPeter Brune @*/
1644422a814eSBarry Smith PetscErrorCode  SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result)
16456a388c36SPeter Brune {
16466a388c36SPeter Brune   PetscFunctionBegin;
16475fd66863SKarl Rupp   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1648422a814eSBarry Smith   linesearch->result = result;
16496a388c36SPeter Brune   PetscFunctionReturn(0);
16506a388c36SPeter Brune }
16516a388c36SPeter Brune 
16529bd66eb0SPeter Brune /*@C
1653f1c6b773SPeter Brune    SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.
16549bd66eb0SPeter Brune 
16559bd66eb0SPeter Brune    Input Parameters:
16569bd66eb0SPeter Brune +  snes - nonlinear context obtained from SNESCreate()
16579bd66eb0SPeter Brune .  projectfunc - function for projecting the function to the bounds
16589bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
16599bd66eb0SPeter Brune 
16609bd66eb0SPeter Brune    Logically Collective on SNES
16619bd66eb0SPeter Brune 
16629bd66eb0SPeter Brune    Calling sequence of projectfunc:
16639bd66eb0SPeter Brune .vb
16649bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X)
16659bd66eb0SPeter Brune .ve
16669bd66eb0SPeter Brune 
16679bd66eb0SPeter Brune     Input parameters for projectfunc:
16689bd66eb0SPeter Brune +   snes - nonlinear context
16699bd66eb0SPeter Brune -   X - current solution
16709bd66eb0SPeter Brune 
1671cd7522eaSPeter Brune     Output parameters for projectfunc:
16729bd66eb0SPeter Brune .   X - Projected solution
16739bd66eb0SPeter Brune 
16749bd66eb0SPeter Brune    Calling sequence of normfunc:
16759bd66eb0SPeter Brune .vb
16769bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm)
16779bd66eb0SPeter Brune .ve
16789bd66eb0SPeter Brune 
1679cd7522eaSPeter Brune     Input parameters for normfunc:
16809bd66eb0SPeter Brune +   snes - nonlinear context
16819bd66eb0SPeter Brune .   X - current solution
16829bd66eb0SPeter Brune -   F - current residual
16839bd66eb0SPeter Brune 
1684cd7522eaSPeter Brune     Output parameters for normfunc:
16859bd66eb0SPeter Brune .   fnorm - VI-specific norm of the function
16869bd66eb0SPeter Brune 
1687cd7522eaSPeter Brune     Notes:
1688cd7522eaSPeter Brune     The VI solvers require projection of the solution to the feasible set.  projectfunc should implement this.
1689cd7522eaSPeter Brune 
1690cd7522eaSPeter Brune     The VI solvers require special evaluation of the function norm such that the norm is only calculated
1691cd7522eaSPeter Brune     on the inactive set.  This should be implemented by normfunc.
16929bd66eb0SPeter Brune 
16939bd66eb0SPeter Brune     Level: developer
16949bd66eb0SPeter Brune 
1695f1c6b773SPeter Brune .seealso: SNESLineSearchGetVIFunctions(), SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck()
16969bd66eb0SPeter Brune @*/
169725acbd8eSLisandro Dalcin PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc)
16989bd66eb0SPeter Brune {
16999bd66eb0SPeter Brune   PetscFunctionBegin;
1700f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
17019bd66eb0SPeter Brune   if (projectfunc) linesearch->ops->viproject = projectfunc;
17029bd66eb0SPeter Brune   if (normfunc) linesearch->ops->vinorm = normfunc;
17039bd66eb0SPeter Brune   PetscFunctionReturn(0);
17049bd66eb0SPeter Brune }
17059bd66eb0SPeter Brune 
17069bd66eb0SPeter Brune /*@C
1707f1c6b773SPeter Brune    SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.
17089bd66eb0SPeter Brune 
17099bd66eb0SPeter Brune    Input Parameters:
1710907376e6SBarry Smith .  linesearch - the line search context, obtain with SNESGetLineSearch()
17119bd66eb0SPeter Brune 
17129bd66eb0SPeter Brune    Output Parameters:
17139bd66eb0SPeter Brune +  projectfunc - function for projecting the function to the bounds
17149bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
17159bd66eb0SPeter Brune 
17169bd66eb0SPeter Brune    Logically Collective on SNES
17179bd66eb0SPeter Brune 
17189bd66eb0SPeter Brune     Level: developer
17199bd66eb0SPeter Brune 
1720f1c6b773SPeter Brune .seealso: SNESLineSearchSetVIFunctions(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
17219bd66eb0SPeter Brune @*/
172225acbd8eSLisandro Dalcin PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc)
17239bd66eb0SPeter Brune {
17249bd66eb0SPeter Brune   PetscFunctionBegin;
17259bd66eb0SPeter Brune   if (projectfunc) *projectfunc = linesearch->ops->viproject;
17269bd66eb0SPeter Brune   if (normfunc) *normfunc = linesearch->ops->vinorm;
17279bd66eb0SPeter Brune   PetscFunctionReturn(0);
17289bd66eb0SPeter Brune }
17299bd66eb0SPeter Brune 
1730bf7f4e0aSPeter Brune /*@C
17311c84c290SBarry Smith   SNESLineSearchRegister - See SNESLineSearchRegister()
1732bf7f4e0aSPeter Brune 
1733bf7f4e0aSPeter Brune   Level: advanced
1734bf7f4e0aSPeter Brune @*/
1735bdf89e91SBarry Smith PetscErrorCode  SNESLineSearchRegister(const char sname[],PetscErrorCode (*function)(SNESLineSearch))
1736bf7f4e0aSPeter Brune {
1737bf7f4e0aSPeter Brune   PetscErrorCode ierr;
1738bf7f4e0aSPeter Brune 
1739bf7f4e0aSPeter Brune   PetscFunctionBegin;
17401d36bdfdSBarry Smith   ierr = SNESInitializePackage();CHKERRQ(ierr);
1741a240a19fSJed Brown   ierr = PetscFunctionListAdd(&SNESLineSearchList,sname,function);CHKERRQ(ierr);
1742bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1743bf7f4e0aSPeter Brune }
1744