xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision 1f8196a2a983893c19b0b7a69c4a23cd5726c506)
1af0996ceSBarry Smith #include <petsc/private/linesearchimpl.h> /*I "petscsnes.h" I*/
2bf7f4e0aSPeter Brune 
3f1c6b773SPeter Brune PetscBool         SNESLineSearchRegisterAllCalled = PETSC_FALSE;
40298fd71SBarry Smith PetscFunctionList SNESLineSearchList              = NULL;
5bf7f4e0aSPeter Brune 
6f1c6b773SPeter Brune PetscClassId  SNESLINESEARCH_CLASSID;
757a83faaSBarry Smith PetscLogEvent SNESLINESEARCH_Apply;
8bf7f4e0aSPeter Brune 
9dcf2fd19SBarry Smith /*@
10dcf2fd19SBarry Smith    SNESLineSearchMonitorCancel - Clears all the monitor functions for a SNESLineSearch object.
11dcf2fd19SBarry Smith 
12dcf2fd19SBarry Smith    Logically Collective on SNESLineSearch
13dcf2fd19SBarry Smith 
14dcf2fd19SBarry Smith    Input Parameters:
15dcf2fd19SBarry Smith .  ls - the SNESLineSearch context
16dcf2fd19SBarry Smith 
17dcf2fd19SBarry Smith    Options Database Key:
18dcf2fd19SBarry Smith .  -snes_linesearch_monitor_cancel - cancels all monitors that have been hardwired
19dcf2fd19SBarry Smith     into a code by calls to SNESLineSearchMonitorSet(), but does not cancel those
20dcf2fd19SBarry Smith     set via the options database
21dcf2fd19SBarry Smith 
22dcf2fd19SBarry Smith    Notes:
23dcf2fd19SBarry Smith    There is no way to clear one specific monitor from a SNESLineSearch object.
24dcf2fd19SBarry Smith 
25dcf2fd19SBarry Smith    This does not clear the monitor set with SNESLineSearchSetDefaultMonitor() use SNESLineSearchSetDefaultMonitor(ls,NULL) to cancel
26dcf2fd19SBarry Smith    that one.
27dcf2fd19SBarry Smith 
28dcf2fd19SBarry Smith    Level: intermediate
29dcf2fd19SBarry Smith 
30dcf2fd19SBarry Smith .keywords: SNESLineSearch, nonlinear, set, monitor
31dcf2fd19SBarry Smith 
32dcf2fd19SBarry Smith .seealso: SNESLineSearchMonitorDefault(), SNESLineSearchMonitorSet()
33dcf2fd19SBarry Smith @*/
34dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchMonitorCancel(SNESLineSearch ls)
35dcf2fd19SBarry Smith {
36dcf2fd19SBarry Smith   PetscErrorCode ierr;
37dcf2fd19SBarry Smith   PetscInt       i;
38dcf2fd19SBarry Smith 
39dcf2fd19SBarry Smith   PetscFunctionBegin;
40dcf2fd19SBarry Smith   PetscValidHeaderSpecific(ls,SNESLINESEARCH_CLASSID,1);
41dcf2fd19SBarry Smith   for (i=0; i<ls->numbermonitors; i++) {
42dcf2fd19SBarry Smith     if (ls->monitordestroy[i]) {
43dcf2fd19SBarry Smith       ierr = (*ls->monitordestroy[i])(&ls->monitorcontext[i]);CHKERRQ(ierr);
44dcf2fd19SBarry Smith     }
45dcf2fd19SBarry Smith   }
46dcf2fd19SBarry Smith   ls->numbermonitors = 0;
47dcf2fd19SBarry Smith   PetscFunctionReturn(0);
48dcf2fd19SBarry Smith }
49dcf2fd19SBarry Smith 
50dcf2fd19SBarry Smith /*@
51dcf2fd19SBarry Smith    SNESLineSearchMonitor - runs the user provided monitor routines, if they exist
52dcf2fd19SBarry Smith 
53dcf2fd19SBarry Smith    Collective on SNES
54dcf2fd19SBarry Smith 
55dcf2fd19SBarry Smith    Input Parameters:
56dcf2fd19SBarry Smith .  ls - the linesearch object
57dcf2fd19SBarry Smith 
58dcf2fd19SBarry Smith    Notes:
59dcf2fd19SBarry Smith    This routine is called by the SNES implementations.
60dcf2fd19SBarry Smith    It does not typically need to be called by the user.
61dcf2fd19SBarry Smith 
62dcf2fd19SBarry Smith    Level: developer
63dcf2fd19SBarry Smith 
64dcf2fd19SBarry Smith .seealso: SNESLineSearchMonitorSet()
65dcf2fd19SBarry Smith @*/
66dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchMonitor(SNESLineSearch ls)
67dcf2fd19SBarry Smith {
68dcf2fd19SBarry Smith   PetscErrorCode ierr;
69dcf2fd19SBarry Smith   PetscInt       i,n = ls->numbermonitors;
70dcf2fd19SBarry Smith 
71dcf2fd19SBarry Smith   PetscFunctionBegin;
72dcf2fd19SBarry Smith   for (i=0; i<n; i++) {
73dcf2fd19SBarry Smith     ierr = (*ls->monitorftns[i])(ls,ls->monitorcontext[i]);CHKERRQ(ierr);
74dcf2fd19SBarry Smith   }
75dcf2fd19SBarry Smith   PetscFunctionReturn(0);
76dcf2fd19SBarry Smith }
77dcf2fd19SBarry Smith 
78dcf2fd19SBarry Smith /*@C
79dcf2fd19SBarry Smith    SNESLineSearchMonitorSet - Sets an ADDITIONAL function that is to be used at every
80dcf2fd19SBarry Smith    iteration of the nonlinear solver to display the iteration's
81dcf2fd19SBarry Smith    progress.
82dcf2fd19SBarry Smith 
83dcf2fd19SBarry Smith    Logically Collective on SNESLineSearch
84dcf2fd19SBarry Smith 
85dcf2fd19SBarry Smith    Input Parameters:
86dcf2fd19SBarry Smith +  ls - the SNESLineSearch context
87dcf2fd19SBarry Smith .  f - the monitor function
88dcf2fd19SBarry Smith .  mctx - [optional] user-defined context for private data for the
89dcf2fd19SBarry Smith           monitor routine (use NULL if no context is desired)
90dcf2fd19SBarry Smith -  monitordestroy - [optional] routine that frees monitor context
91dcf2fd19SBarry Smith           (may be NULL)
92dcf2fd19SBarry Smith 
93dcf2fd19SBarry Smith    Notes:
94dcf2fd19SBarry Smith    Several different monitoring routines may be set by calling
95dcf2fd19SBarry Smith    SNESLineSearchMonitorSet() multiple times; all will be called in the
96dcf2fd19SBarry Smith    order in which they were set.
97dcf2fd19SBarry Smith 
98dcf2fd19SBarry Smith    Fortran notes: Only a single monitor function can be set for each SNESLineSearch object
99dcf2fd19SBarry Smith 
100dcf2fd19SBarry Smith    Level: intermediate
101dcf2fd19SBarry Smith 
102dcf2fd19SBarry Smith .keywords: SNESLineSearch, nonlinear, set, monitor
103dcf2fd19SBarry Smith 
104dcf2fd19SBarry Smith .seealso: SNESLineSearchMonitorDefault(), SNESLineSearchMonitorCancel()
105dcf2fd19SBarry Smith @*/
106dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchMonitorSet(SNESLineSearch ls,PetscErrorCode (*f)(SNESLineSearch,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
107dcf2fd19SBarry Smith {
10878064530SBarry Smith   PetscErrorCode ierr;
10978064530SBarry Smith   PetscInt       i;
11078064530SBarry Smith   PetscBool      identical;
11178064530SBarry Smith 
112dcf2fd19SBarry Smith   PetscFunctionBegin;
113dcf2fd19SBarry Smith   PetscValidHeaderSpecific(ls,SNESLINESEARCH_CLASSID,1);
11478064530SBarry Smith   for (i=0; i<ls->numbermonitors;i++) {
11578064530SBarry Smith     ierr = PetscMonitorCompare((PetscErrorCode (*)(void))f,mctx,monitordestroy,(PetscErrorCode (*)(void))ls->monitorftns[i],ls->monitorcontext[i],ls->monitordestroy[i],&identical);CHKERRQ(ierr);
11678064530SBarry Smith     if (identical) PetscFunctionReturn(0);
11778064530SBarry Smith   }
118dcf2fd19SBarry Smith   if (ls->numbermonitors >= MAXSNESLSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
119dcf2fd19SBarry Smith   ls->monitorftns[ls->numbermonitors]          = f;
120dcf2fd19SBarry Smith   ls->monitordestroy[ls->numbermonitors]   = monitordestroy;
121dcf2fd19SBarry Smith   ls->monitorcontext[ls->numbermonitors++] = (void*)mctx;
122dcf2fd19SBarry Smith   PetscFunctionReturn(0);
123dcf2fd19SBarry Smith }
124dcf2fd19SBarry Smith 
125dcf2fd19SBarry Smith /*@C
126dcf2fd19SBarry Smith    SNESLineSearchMonitorSolutionUpdate - Monitors each update a new function value the linesearch tries
127dcf2fd19SBarry Smith 
128dcf2fd19SBarry Smith    Collective on SNESLineSearch
129dcf2fd19SBarry Smith 
130dcf2fd19SBarry Smith    Input Parameters:
131dcf2fd19SBarry Smith +  ls - the SNES linesearch object
132d12e167eSBarry Smith -  vf - the context for the monitor, in this case it is an ASCII PetscViewer and format
133dcf2fd19SBarry Smith 
134dcf2fd19SBarry Smith    Level: intermediate
135dcf2fd19SBarry Smith 
136dcf2fd19SBarry Smith .keywords: SNES, nonlinear, default, monitor, norm
137dcf2fd19SBarry Smith 
138dcf2fd19SBarry Smith .seealso: SNESMonitorSet(), SNESMonitorSolution()
139dcf2fd19SBarry Smith @*/
140d12e167eSBarry Smith PetscErrorCode  SNESLineSearchMonitorSolutionUpdate(SNESLineSearch ls,PetscViewerAndFormat *vf)
141dcf2fd19SBarry Smith {
142dcf2fd19SBarry Smith   PetscErrorCode ierr;
143d12e167eSBarry Smith   PetscViewer    viewer = vf->viewer;
144dcf2fd19SBarry Smith   Vec            Y,W,G;
145dcf2fd19SBarry Smith 
146dcf2fd19SBarry Smith   PetscFunctionBegin;
147dcf2fd19SBarry Smith   ierr = SNESLineSearchGetVecs(ls,NULL,NULL,&Y,&W,&G);CHKERRQ(ierr);
148d12e167eSBarry Smith   ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr);
149dcf2fd19SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"LineSearch attempted update to solution \n");CHKERRQ(ierr);
150dcf2fd19SBarry Smith   ierr = VecView(Y,viewer);CHKERRQ(ierr);
151dcf2fd19SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"LineSearch attempted new solution \n");CHKERRQ(ierr);
152dcf2fd19SBarry Smith   ierr = VecView(W,viewer);CHKERRQ(ierr);
153dcf2fd19SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"LineSearch attempted updated function value\n");CHKERRQ(ierr);
154dcf2fd19SBarry Smith   ierr = VecView(G,viewer);CHKERRQ(ierr);
155d12e167eSBarry Smith   ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
156dcf2fd19SBarry Smith   PetscFunctionReturn(0);
157dcf2fd19SBarry Smith }
158dcf2fd19SBarry Smith 
159f40b411bSPeter Brune /*@
160cd7522eaSPeter Brune    SNESLineSearchCreate - Creates the line search context.
161f40b411bSPeter Brune 
162cd7522eaSPeter Brune    Logically Collective on Comm
163f40b411bSPeter Brune 
164f40b411bSPeter Brune    Input Parameters:
165cd7522eaSPeter Brune .  comm - MPI communicator for the line search (typically from the associated SNES context).
166f40b411bSPeter Brune 
167f40b411bSPeter Brune    Output Parameters:
1688e557f58SPeter Brune .  outlinesearch - the new linesearch context
169f40b411bSPeter Brune 
170162e0bf5SPeter Brune    Level: developer
171162e0bf5SPeter Brune 
172162e0bf5SPeter Brune    Notes:
173162e0bf5SPeter Brune    The preferred calling sequence for users is to use SNESGetLineSearch() to acquire the SNESLineSearch instance
174162e0bf5SPeter Brune    already associated with the SNES.  This function is for developer use.
175f40b411bSPeter Brune 
176cd7522eaSPeter Brune .keywords: LineSearch, create, context
177f40b411bSPeter Brune 
178162e0bf5SPeter Brune .seealso: LineSearchDestroy(), SNESGetLineSearch()
179f40b411bSPeter Brune @*/
180f40b411bSPeter Brune 
181bf388a1fSBarry Smith PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch)
182bf388a1fSBarry Smith {
183bf7f4e0aSPeter Brune   PetscErrorCode ierr;
184f1c6b773SPeter Brune   SNESLineSearch linesearch;
185bf388a1fSBarry Smith 
186bf7f4e0aSPeter Brune   PetscFunctionBegin;
187ea5d4fccSPeter Brune   PetscValidPointer(outlinesearch,2);
188f34a81feSMatthew G. Knepley   ierr = SNESInitializePackage();CHKERRQ(ierr);
1890298fd71SBarry Smith   *outlinesearch = NULL;
190f5af7f23SKarl Rupp 
19173107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(linesearch,SNESLINESEARCH_CLASSID, "SNESLineSearch","Linesearch","SNESLineSearch",comm,SNESLineSearchDestroy,SNESLineSearchView);CHKERRQ(ierr);
192bf7f4e0aSPeter Brune 
1930298fd71SBarry Smith   linesearch->vec_sol_new  = NULL;
1940298fd71SBarry Smith   linesearch->vec_func_new = NULL;
1950298fd71SBarry Smith   linesearch->vec_sol      = NULL;
1960298fd71SBarry Smith   linesearch->vec_func     = NULL;
1970298fd71SBarry Smith   linesearch->vec_update   = NULL;
1989bd66eb0SPeter Brune 
199bf7f4e0aSPeter Brune   linesearch->lambda       = 1.0;
200bf7f4e0aSPeter Brune   linesearch->fnorm        = 1.0;
201bf7f4e0aSPeter Brune   linesearch->ynorm        = 1.0;
202bf7f4e0aSPeter Brune   linesearch->xnorm        = 1.0;
203422a814eSBarry Smith   linesearch->result       = SNES_LINESEARCH_SUCCEEDED;
204bf7f4e0aSPeter Brune   linesearch->norms        = PETSC_TRUE;
205bf7f4e0aSPeter Brune   linesearch->keeplambda   = PETSC_FALSE;
206bf7f4e0aSPeter Brune   linesearch->damping      = 1.0;
207bf7f4e0aSPeter Brune   linesearch->maxstep      = 1e8;
208bf7f4e0aSPeter Brune   linesearch->steptol      = 1e-12;
209516fe3c3SPeter Brune   linesearch->rtol         = 1e-8;
210516fe3c3SPeter Brune   linesearch->atol         = 1e-15;
211516fe3c3SPeter Brune   linesearch->ltol         = 1e-8;
2120298fd71SBarry Smith   linesearch->precheckctx  = NULL;
2130298fd71SBarry Smith   linesearch->postcheckctx = NULL;
21459405d5eSPeter Brune   linesearch->max_its      = 1;
215bf7f4e0aSPeter Brune   linesearch->setupcalled  = PETSC_FALSE;
216bf7f4e0aSPeter Brune   *outlinesearch           = linesearch;
217bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
218bf7f4e0aSPeter Brune }
219bf7f4e0aSPeter Brune 
220f40b411bSPeter Brune /*@
22178bcb3b5SPeter Brune    SNESLineSearchSetUp - Prepares the line search for being applied by allocating
22278bcb3b5SPeter Brune    any required vectors.
223f40b411bSPeter Brune 
224cd7522eaSPeter Brune    Collective on SNESLineSearch
225f40b411bSPeter Brune 
226f40b411bSPeter Brune    Input Parameters:
227f40b411bSPeter Brune .  linesearch - The LineSearch instance.
228f40b411bSPeter Brune 
229cd7522eaSPeter Brune    Notes:
230f190f2fcSBarry Smith    For most cases, this needn't be called by users or outside of SNESLineSearchApply().
231cd7522eaSPeter Brune    The only current case where this is called outside of this is for the VI
23278bcb3b5SPeter Brune    solvers, which modify the solution and work vectors before the first call
233cd7522eaSPeter Brune    of SNESLineSearchApply, requiring the SNESLineSearch work vectors to be
234cd7522eaSPeter Brune    allocated upfront.
235cd7522eaSPeter Brune 
23678bcb3b5SPeter Brune    Level: advanced
237f40b411bSPeter Brune 
238f1c6b773SPeter Brune .keywords: SNESLineSearch, SetUp
239f40b411bSPeter Brune 
240f1c6b773SPeter Brune .seealso: SNESLineSearchReset()
241f40b411bSPeter Brune @*/
242f40b411bSPeter Brune 
243bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch)
244bf388a1fSBarry Smith {
245bf7f4e0aSPeter Brune   PetscErrorCode ierr;
246bf388a1fSBarry Smith 
247bf7f4e0aSPeter Brune   PetscFunctionBegin;
248bf7f4e0aSPeter Brune   if (!((PetscObject)linesearch)->type_name) {
2491a4f838cSPeter Brune     ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr);
250bf7f4e0aSPeter Brune   }
251bf7f4e0aSPeter Brune   if (!linesearch->setupcalled) {
2529bd66eb0SPeter Brune     if (!linesearch->vec_sol_new) {
253bf7f4e0aSPeter Brune       ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new);CHKERRQ(ierr);
2549bd66eb0SPeter Brune     }
2559bd66eb0SPeter Brune     if (!linesearch->vec_func_new) {
2569bd66eb0SPeter Brune       ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new);CHKERRQ(ierr);
2579bd66eb0SPeter Brune     }
258bf7f4e0aSPeter Brune     if (linesearch->ops->setup) {
259bf7f4e0aSPeter Brune       ierr = (*linesearch->ops->setup)(linesearch);CHKERRQ(ierr);
260bf7f4e0aSPeter Brune     }
261ed07d7d7SPeter Brune     if (!linesearch->ops->snesfunc) {ierr = SNESLineSearchSetFunction(linesearch,SNESComputeFunction);CHKERRQ(ierr);}
262bf7f4e0aSPeter Brune     linesearch->lambda      = linesearch->damping;
263bf7f4e0aSPeter Brune     linesearch->setupcalled = PETSC_TRUE;
264bf7f4e0aSPeter Brune   }
265bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
266bf7f4e0aSPeter Brune }
267bf7f4e0aSPeter Brune 
268f40b411bSPeter Brune 
269f40b411bSPeter Brune /*@
270f190f2fcSBarry Smith    SNESLineSearchReset - Undoes the SNESLineSearchSetUp() and deletes any Vecs or Mats allocated by the line search.
271f40b411bSPeter Brune 
272f1c6b773SPeter Brune    Collective on SNESLineSearch
273f40b411bSPeter Brune 
274f40b411bSPeter Brune    Input Parameters:
275f40b411bSPeter Brune .  linesearch - The LineSearch instance.
276f40b411bSPeter Brune 
277f190f2fcSBarry Smith    Notes: Usually only called by SNESReset()
278f190f2fcSBarry Smith 
279f190f2fcSBarry Smith    Level: developer
280f40b411bSPeter Brune 
281cd7522eaSPeter Brune .keywords: SNESLineSearch, Reset
282f40b411bSPeter Brune 
283f1c6b773SPeter Brune .seealso: SNESLineSearchSetUp()
284f40b411bSPeter Brune @*/
285f40b411bSPeter Brune 
286bf388a1fSBarry Smith PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch)
287bf388a1fSBarry Smith {
288bf7f4e0aSPeter Brune   PetscErrorCode ierr;
289bf388a1fSBarry Smith 
290bf7f4e0aSPeter Brune   PetscFunctionBegin;
291f5af7f23SKarl Rupp   if (linesearch->ops->reset) (*linesearch->ops->reset)(linesearch);
292f5af7f23SKarl Rupp 
293bf7f4e0aSPeter Brune   ierr = VecDestroy(&linesearch->vec_sol_new);CHKERRQ(ierr);
294bf7f4e0aSPeter Brune   ierr = VecDestroy(&linesearch->vec_func_new);CHKERRQ(ierr);
295bf7f4e0aSPeter Brune 
296bf7f4e0aSPeter Brune   ierr = VecDestroyVecs(linesearch->nwork, &linesearch->work);CHKERRQ(ierr);
297f5af7f23SKarl Rupp 
298bf7f4e0aSPeter Brune   linesearch->nwork       = 0;
299bf7f4e0aSPeter Brune   linesearch->setupcalled = PETSC_FALSE;
300bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
301bf7f4e0aSPeter Brune }
302bf7f4e0aSPeter Brune 
303ed07d7d7SPeter Brune /*@C
304f190f2fcSBarry Smith    SNESLineSearchSetFunction - Sets the function evaluation used by the SNES line search
305ed07d7d7SPeter Brune 
306ed07d7d7SPeter Brune    Input Parameters:
307ed07d7d7SPeter Brune .  linesearch - the SNESLineSearch context
308f190f2fcSBarry Smith +  func       - function evaluation routine
309ed07d7d7SPeter Brune 
310ed07d7d7SPeter Brune    Level: developer
311ed07d7d7SPeter Brune 
312f190f2fcSBarry Smith    Notes: This is used internally by PETSc and not called by users
313f190f2fcSBarry Smith 
314ed07d7d7SPeter Brune .keywords: get, linesearch, pre-check
315ed07d7d7SPeter Brune 
316f190f2fcSBarry Smith .seealso: SNESSetFunction()
317ed07d7d7SPeter Brune @*/
318ed07d7d7SPeter Brune PetscErrorCode  SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES,Vec,Vec))
319ed07d7d7SPeter Brune {
320ed07d7d7SPeter Brune   PetscFunctionBegin;
321ed07d7d7SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
322ed07d7d7SPeter Brune   linesearch->ops->snesfunc = func;
323ed07d7d7SPeter Brune   PetscFunctionReturn(0);
324ed07d7d7SPeter Brune }
325ed07d7d7SPeter Brune 
326ed07d7d7SPeter Brune 
3276b2b7091SBarry Smith /*MC
328f190f2fcSBarry Smith     SNESLineSearchPreCheckFunction - form of function passed to check the search direction before line search is called
3296b2b7091SBarry Smith 
3306b2b7091SBarry Smith      Synopsis:
331aaa7dc30SBarry Smith      #include <petscsnes.h>
3326b2b7091SBarry Smith      SNESLineSearchPreCheckFunction(SNESLineSearch snes,Vec x,Vec y, PetscBool *changed);
3336b2b7091SBarry Smith 
3346b2b7091SBarry Smith        Input Parameters:
3356b2b7091SBarry Smith +      x - solution vector
3366b2b7091SBarry Smith .      y - search direction vector
3376b2b7091SBarry Smith -      changed - flag to indicate the precheck changed x or y.
3386b2b7091SBarry Smith 
339f190f2fcSBarry Smith      Note: This is NOTE a PETSc function, rather it documents the calling sequence of functions passed to SNESLineSearchSetPreCheck()
340f190f2fcSBarry Smith            and SNESLineSearchGetPreCheck()
341f190f2fcSBarry Smith 
342878cb397SSatish Balay    Level: advanced
343878cb397SSatish Balay 
344f190f2fcSBarry Smith .seealso:   SNESLineSearchSetPreCheck(), SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck()
3456b2b7091SBarry Smith M*/
3466b2b7091SBarry Smith 
34786d74e61SPeter Brune /*@C
348f190f2fcSBarry Smith    SNESLineSearchSetPreCheck - Sets a user function that is called after the initial search direction has been computed but
349df3898eeSBarry Smith          before the line search routine has been applied. Allows the user to adjust the result of (usually a linear solve) that
350f190f2fcSBarry Smith          determined the search direction.
35186d74e61SPeter Brune 
352f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
35386d74e61SPeter Brune 
35486d74e61SPeter Brune    Input Parameters:
355f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
356f190f2fcSBarry Smith .  func - [optional] function evaluation routine, see SNESLineSearchPreCheckFunction for the calling sequence
357f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
35886d74e61SPeter Brune 
35986d74e61SPeter Brune    Level: intermediate
36086d74e61SPeter Brune 
36186d74e61SPeter Brune .keywords: set, linesearch, pre-check
36286d74e61SPeter Brune 
363f190f2fcSBarry Smith .seealso: SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
36486d74e61SPeter Brune @*/
365f190f2fcSBarry Smith PetscErrorCode  SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void *ctx)
36686d74e61SPeter Brune {
3679bd66eb0SPeter Brune   PetscFunctionBegin;
368f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
369f190f2fcSBarry Smith   if (func) linesearch->ops->precheck = func;
37086d74e61SPeter Brune   if (ctx) linesearch->precheckctx = ctx;
37186d74e61SPeter Brune   PetscFunctionReturn(0);
37286d74e61SPeter Brune }
37386d74e61SPeter Brune 
37486d74e61SPeter Brune /*@C
37553deb899SBarry Smith    SNESLineSearchGetPreCheck - Gets the pre-check function for the line search routine.
37686d74e61SPeter Brune 
37786d74e61SPeter Brune    Input Parameters:
378f1c6b773SPeter Brune .  linesearch - the SNESLineSearch context
37986d74e61SPeter Brune 
38086d74e61SPeter Brune    Output Parameters:
381f190f2fcSBarry Smith +  func       - [optional] function evaluation routine, see SNESLineSearchPreCheckFunction for calling sequence
382f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
38386d74e61SPeter Brune 
38486d74e61SPeter Brune    Level: intermediate
38586d74e61SPeter Brune 
38686d74e61SPeter Brune .keywords: get, linesearch, pre-check
38786d74e61SPeter Brune 
388f1c6b773SPeter Brune .seealso: SNESLineSearchGetPostCheck(), SNESLineSearchSetPreCheck()
38986d74e61SPeter Brune @*/
390f190f2fcSBarry Smith PetscErrorCode  SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void **ctx)
39186d74e61SPeter Brune {
3929bd66eb0SPeter Brune   PetscFunctionBegin;
393f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
394f190f2fcSBarry Smith   if (func) *func = linesearch->ops->precheck;
39586d74e61SPeter Brune   if (ctx) *ctx = linesearch->precheckctx;
39686d74e61SPeter Brune   PetscFunctionReturn(0);
39786d74e61SPeter Brune }
39886d74e61SPeter Brune 
3996b2b7091SBarry Smith /*MC
400f190f2fcSBarry Smith     SNESLineSearchPostCheckFunction - form of function that is called after line search is complete
4016b2b7091SBarry Smith 
4026b2b7091SBarry Smith      Synopsis:
403aaa7dc30SBarry Smith      #include <petscsnes.h>
4046b2b7091SBarry Smith      SNESLineSearchPostheckFunction(SNESLineSearch linesearch,Vec x,Vec y,  Vec w, *changed_y, PetscBool *changed_w);
4056b2b7091SBarry Smith 
4066b2b7091SBarry Smith      Input Parameters:
4076b2b7091SBarry Smith +      x - old solution vector
4086b2b7091SBarry Smith .      y - search direction vector
4096b2b7091SBarry Smith .      w - new solution vector
4106b2b7091SBarry Smith .      changed_y - indicates that the line search changed y
4116b2b7091SBarry Smith -      changed_w - indicates that the line search changed w
4126b2b7091SBarry Smith 
413f190f2fcSBarry Smith      Note: This is NOTE a PETSc function, rather it documents the calling sequence of functions passed to SNESLineSearchSetPostCheck()
414f190f2fcSBarry Smith            and SNESLineSearchGetPostCheck()
415f190f2fcSBarry Smith 
416878cb397SSatish Balay    Level: advanced
4176b2b7091SBarry Smith 
418f190f2fcSBarry Smith .seealso:   SNESLineSearchSetPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPreCheck(), SNESLineSearchGetPostCheck()
4196b2b7091SBarry Smith M*/
4206b2b7091SBarry Smith 
42186d74e61SPeter Brune /*@C
422f190f2fcSBarry Smith    SNESLineSearchSetPostCheck - Sets a user function that is called after the line search has been applied to determine the step
423f190f2fcSBarry Smith        direction and length. Allows the user a chance to change or override the decision of the line search routine
42486d74e61SPeter Brune 
425f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
42686d74e61SPeter Brune 
42786d74e61SPeter Brune    Input Parameters:
428f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
429f190f2fcSBarry Smith .  func - [optional] function evaluation routine, see SNESLineSearchPostCheckFunction for the calling sequence
430f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
43186d74e61SPeter Brune 
43286d74e61SPeter Brune    Level: intermediate
43386d74e61SPeter Brune 
43486d74e61SPeter Brune .keywords: set, linesearch, post-check
43586d74e61SPeter Brune 
436f1c6b773SPeter Brune .seealso: SNESLineSearchSetPreCheck()
43786d74e61SPeter Brune @*/
438f190f2fcSBarry Smith PetscErrorCode  SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void *ctx)
43986d74e61SPeter Brune {
44086d74e61SPeter Brune   PetscFunctionBegin;
441f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
442f190f2fcSBarry Smith   if (func) linesearch->ops->postcheck = func;
44386d74e61SPeter Brune   if (ctx) linesearch->postcheckctx = ctx;
44486d74e61SPeter Brune   PetscFunctionReturn(0);
44586d74e61SPeter Brune }
44686d74e61SPeter Brune 
44786d74e61SPeter Brune /*@C
448f1c6b773SPeter Brune    SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine.
44986d74e61SPeter Brune 
45086d74e61SPeter Brune    Input Parameters:
451f1c6b773SPeter Brune .  linesearch - the SNESLineSearch context
45286d74e61SPeter Brune 
45386d74e61SPeter Brune    Output Parameters:
454f190f2fcSBarry Smith +  func - [optional] function evaluation routine, see for the calling sequence SNESLineSearchPostCheckFunction
455f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
45686d74e61SPeter Brune 
45786d74e61SPeter Brune    Level: intermediate
45886d74e61SPeter Brune 
45986d74e61SPeter Brune .keywords: get, linesearch, post-check
46086d74e61SPeter Brune 
461f1c6b773SPeter Brune .seealso: SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck()
46286d74e61SPeter Brune @*/
463f190f2fcSBarry Smith PetscErrorCode  SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void **ctx)
46486d74e61SPeter Brune {
4659bd66eb0SPeter Brune   PetscFunctionBegin;
466f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
467f190f2fcSBarry Smith   if (func) *func = linesearch->ops->postcheck;
46886d74e61SPeter Brune   if (ctx) *ctx = linesearch->postcheckctx;
46986d74e61SPeter Brune   PetscFunctionReturn(0);
47086d74e61SPeter Brune }
47186d74e61SPeter Brune 
472f40b411bSPeter Brune /*@
473f1c6b773SPeter Brune    SNESLineSearchPreCheck - Prepares the line search for being applied.
474f40b411bSPeter Brune 
475cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
476f40b411bSPeter Brune 
477f40b411bSPeter Brune    Input Parameters:
4787b1df9c1SPeter Brune +  linesearch - The linesearch instance.
4797b1df9c1SPeter Brune .  X - The current solution
4807b1df9c1SPeter Brune -  Y - The step direction
481f40b411bSPeter Brune 
482f40b411bSPeter Brune    Output Parameters:
4838e557f58SPeter Brune .  changed - Indicator that the precheck routine has changed anything
484f40b411bSPeter Brune 
485f190f2fcSBarry Smith    Level: developer
486f40b411bSPeter Brune 
487f1c6b773SPeter Brune .keywords: SNESLineSearch, Create
488f40b411bSPeter Brune 
489f1c6b773SPeter Brune .seealso: SNESLineSearchPostCheck()
490f40b411bSPeter Brune @*/
4917b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed)
492bf7f4e0aSPeter Brune {
493bf7f4e0aSPeter Brune   PetscErrorCode ierr;
4945fd66863SKarl Rupp 
495bf7f4e0aSPeter Brune   PetscFunctionBegin;
496bf7f4e0aSPeter Brune   *changed = PETSC_FALSE;
4976b2b7091SBarry Smith   if (linesearch->ops->precheck) {
4986b2b7091SBarry Smith     ierr = (*linesearch->ops->precheck)(linesearch, X, Y, changed, linesearch->precheckctx);CHKERRQ(ierr);
49938bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed,4);
500bf7f4e0aSPeter Brune   }
501bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
502bf7f4e0aSPeter Brune }
503bf7f4e0aSPeter Brune 
504f40b411bSPeter Brune /*@
505f1c6b773SPeter Brune    SNESLineSearchPostCheck - Prepares the line search for being applied.
506f40b411bSPeter Brune 
507cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
508f40b411bSPeter Brune 
509f40b411bSPeter Brune    Input Parameters:
5107b1df9c1SPeter Brune +  linesearch - The linesearch context
5117b1df9c1SPeter Brune .  X - The last solution
5127b1df9c1SPeter Brune .  Y - The step direction
5137b1df9c1SPeter Brune -  W - The updated solution, W = X + lambda*Y for some lambda
514f40b411bSPeter Brune 
515f40b411bSPeter Brune    Output Parameters:
51678bcb3b5SPeter Brune +  changed_Y - Indicator if the direction Y has been changed.
51778bcb3b5SPeter Brune -  changed_W - Indicator if the new candidate solution W has been changed.
518f40b411bSPeter Brune 
519f190f2fcSBarry Smith    Level: developer
520f40b411bSPeter Brune 
521f1c6b773SPeter Brune .keywords: SNESLineSearch, Create
522f40b411bSPeter Brune 
523f1c6b773SPeter Brune .seealso: SNESLineSearchPreCheck()
524f40b411bSPeter Brune @*/
5257b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch,Vec X,Vec Y,Vec W,PetscBool *changed_Y,PetscBool *changed_W)
526bf7f4e0aSPeter Brune {
527bf7f4e0aSPeter Brune   PetscErrorCode ierr;
528bf388a1fSBarry Smith 
529bf7f4e0aSPeter Brune   PetscFunctionBegin;
530bf7f4e0aSPeter Brune   *changed_Y = PETSC_FALSE;
531bf7f4e0aSPeter Brune   *changed_W = PETSC_FALSE;
5326b2b7091SBarry Smith   if (linesearch->ops->postcheck) {
5336b2b7091SBarry Smith     ierr = (*linesearch->ops->postcheck)(linesearch,X,Y,W,changed_Y,changed_W,linesearch->postcheckctx);CHKERRQ(ierr);
53438bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed_Y,5);
53538bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed_W,6);
53686d74e61SPeter Brune   }
53786d74e61SPeter Brune   PetscFunctionReturn(0);
53886d74e61SPeter Brune }
53986d74e61SPeter Brune 
54086d74e61SPeter Brune /*@C
54186d74e61SPeter Brune    SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration
54286d74e61SPeter Brune 
543cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
54486d74e61SPeter Brune 
54586d74e61SPeter Brune    Input Arguments:
54686d74e61SPeter Brune +  linesearch - linesearch context
54786d74e61SPeter Brune .  X - base state for this step
54886d74e61SPeter Brune .  Y - initial correction
549907376e6SBarry Smith -  ctx - context for this function
55086d74e61SPeter Brune 
55186d74e61SPeter Brune    Output Arguments:
55286d74e61SPeter Brune +  Y - correction, possibly modified
55386d74e61SPeter Brune -  changed - flag indicating that Y was modified
55486d74e61SPeter Brune 
55586d74e61SPeter Brune    Options Database Key:
556cd7522eaSPeter Brune +  -snes_linesearch_precheck_picard - activate this routine
557cd7522eaSPeter Brune -  -snes_linesearch_precheck_picard_angle - angle
55886d74e61SPeter Brune 
55986d74e61SPeter Brune    Level: advanced
56086d74e61SPeter Brune 
56186d74e61SPeter Brune    Notes:
56286d74e61SPeter Brune    This function should be passed to SNESLineSearchSetPreCheck()
56386d74e61SPeter Brune 
56486d74e61SPeter Brune    The justification for this method involves the linear convergence of a Picard iteration
56586d74e61SPeter Brune    so the Picard linearization should be provided in place of the "Jacobian". This correction
56686d74e61SPeter Brune    is generally not useful when using a Newton linearization.
56786d74e61SPeter Brune 
56886d74e61SPeter Brune    Reference:
56986d74e61SPeter Brune    Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology.
57086d74e61SPeter Brune 
57186d74e61SPeter Brune .seealso: SNESLineSearchSetPreCheck()
57286d74e61SPeter Brune @*/
573f1c6b773SPeter Brune PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed,void *ctx)
57486d74e61SPeter Brune {
57586d74e61SPeter Brune   PetscErrorCode ierr;
57686d74e61SPeter Brune   PetscReal      angle = *(PetscReal*)linesearch->precheckctx;
57786d74e61SPeter Brune   Vec            Ylast;
57886d74e61SPeter Brune   PetscScalar    dot;
57986d74e61SPeter Brune   PetscInt       iter;
58086d74e61SPeter Brune   PetscReal      ynorm,ylastnorm,theta,angle_radians;
58186d74e61SPeter Brune   SNES           snes;
58286d74e61SPeter Brune 
58386d74e61SPeter Brune   PetscFunctionBegin;
584f1c6b773SPeter Brune   ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
58586d74e61SPeter Brune   ierr = PetscObjectQuery((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject*)&Ylast);CHKERRQ(ierr);
58686d74e61SPeter Brune   if (!Ylast) {
58786d74e61SPeter Brune     ierr = VecDuplicate(Y,&Ylast);CHKERRQ(ierr);
58886d74e61SPeter Brune     ierr = PetscObjectCompose((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject)Ylast);CHKERRQ(ierr);
58986d74e61SPeter Brune     ierr = PetscObjectDereference((PetscObject)Ylast);CHKERRQ(ierr);
59086d74e61SPeter Brune   }
59186d74e61SPeter Brune   ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr);
59286d74e61SPeter Brune   if (iter < 2) {
59386d74e61SPeter Brune     ierr     = VecCopy(Y,Ylast);CHKERRQ(ierr);
59486d74e61SPeter Brune     *changed = PETSC_FALSE;
59586d74e61SPeter Brune     PetscFunctionReturn(0);
59686d74e61SPeter Brune   }
59786d74e61SPeter Brune 
59886d74e61SPeter Brune   ierr = VecDot(Y,Ylast,&dot);CHKERRQ(ierr);
59986d74e61SPeter Brune   ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr);
60086d74e61SPeter Brune   ierr = VecNorm(Ylast,NORM_2,&ylastnorm);CHKERRQ(ierr);
60186d74e61SPeter Brune   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
602255453a1SBarry Smith   theta         = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0));
60386d74e61SPeter Brune   angle_radians = angle * PETSC_PI / 180.;
60486d74e61SPeter Brune   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
60586d74e61SPeter Brune     /* Modify the step Y */
60686d74e61SPeter Brune     PetscReal alpha,ydiffnorm;
60786d74e61SPeter Brune     ierr  = VecAXPY(Ylast,-1.0,Y);CHKERRQ(ierr);
60886d74e61SPeter Brune     ierr  = VecNorm(Ylast,NORM_2,&ydiffnorm);CHKERRQ(ierr);
60986d74e61SPeter Brune     alpha = ylastnorm / ydiffnorm;
61086d74e61SPeter Brune     ierr  = VecCopy(Y,Ylast);CHKERRQ(ierr);
61186d74e61SPeter Brune     ierr  = VecScale(Y,alpha);CHKERRQ(ierr);
612c69d1a72SBarry 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);
61386d74e61SPeter Brune   } else {
614c69d1a72SBarry 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);
61586d74e61SPeter Brune     ierr     = VecCopy(Y,Ylast);CHKERRQ(ierr);
61686d74e61SPeter Brune     *changed = PETSC_FALSE;
617bf7f4e0aSPeter Brune   }
618bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
619bf7f4e0aSPeter Brune }
620bf7f4e0aSPeter Brune 
621f40b411bSPeter Brune /*@
622cd7522eaSPeter Brune    SNESLineSearchApply - Computes the line-search update.
623f40b411bSPeter Brune 
624f1c6b773SPeter Brune    Collective on SNESLineSearch
625f40b411bSPeter Brune 
626f40b411bSPeter Brune    Input Parameters:
6278e557f58SPeter Brune +  linesearch - The linesearch context
6288e557f58SPeter Brune .  X - The current solution
6298e557f58SPeter Brune .  F - The current function
6308e557f58SPeter Brune .  fnorm - The current norm
6318e557f58SPeter Brune -  Y - The search direction
632f40b411bSPeter Brune 
633f40b411bSPeter Brune    Output Parameters:
6348e557f58SPeter Brune +  X - The new solution
6358e557f58SPeter Brune .  F - The new function
6368e557f58SPeter Brune -  fnorm - The new function norm
637f40b411bSPeter Brune 
638cd7522eaSPeter Brune    Options Database Keys:
639d4c6564cSPatrick Farrell + -snes_linesearch_type - basic, bt, l2, cp, nleqerr, shell
640dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches
6411fe24845SBarry Smith . -snes_linesearch_damping - The linesearch damping parameter, default is 1.0 (no damping)
6421fe24845SBarry Smith . -snes_linesearch_norms   - Turn on/off the linesearch norms computation (SNESLineSearchSetComputeNorms())
6433c7d6663SPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess
6443c7d6663SPeter Brune - -snes_linesearch_max_it - The number of iterations for iterative line searches
645cd7522eaSPeter Brune 
646cd7522eaSPeter Brune    Notes:
647cd7522eaSPeter Brune    This is typically called from within a SNESSolve() implementation in order to
648cd7522eaSPeter Brune    help with convergence of the nonlinear method.  Various SNES types use line searches
649cd7522eaSPeter Brune    in different ways, but the overarching theme is that a line search is used to determine
650cd7522eaSPeter Brune    an optimal damping parameter of a step at each iteration of the method.  Each
651cd7522eaSPeter Brune    application of the line search may invoke SNESComputeFunction several times, and
652cd7522eaSPeter Brune    therefore may be fairly expensive.
653cd7522eaSPeter Brune 
654f40b411bSPeter Brune    Level: Intermediate
655f40b411bSPeter Brune 
656f1c6b773SPeter Brune .keywords: SNESLineSearch, Create
657f40b411bSPeter Brune 
6581fe24845SBarry Smith .seealso: SNESLineSearchCreate(), SNESLineSearchPreCheck(), SNESLineSearchPostCheck(), SNESSolve(), SNESComputeFunction(), SNESLineSearchSetComputeNorms(),
6591fe24845SBarry Smith           SNESLineSearchType, SNESLineSearchSetType()
660f40b411bSPeter Brune @*/
661bf388a1fSBarry Smith PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y)
662bf388a1fSBarry Smith {
663bf7f4e0aSPeter Brune   PetscErrorCode ierr;
664bf7f4e0aSPeter Brune 
665bf388a1fSBarry Smith   PetscFunctionBegin;
666f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
667bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
668bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
669bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
670bf7f4e0aSPeter Brune 
671422a814eSBarry Smith   linesearch->result = SNES_LINESEARCH_SUCCEEDED;
672bf7f4e0aSPeter Brune 
673bf7f4e0aSPeter Brune   linesearch->vec_sol    = X;
674bf7f4e0aSPeter Brune   linesearch->vec_update = Y;
675bf7f4e0aSPeter Brune   linesearch->vec_func   = F;
676bf7f4e0aSPeter Brune 
677f1c6b773SPeter Brune   ierr = SNESLineSearchSetUp(linesearch);CHKERRQ(ierr);
678bf7f4e0aSPeter Brune 
679f5af7f23SKarl Rupp   if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */
680bf7f4e0aSPeter Brune 
681f5af7f23SKarl Rupp   if (fnorm) linesearch->fnorm = *fnorm;
682f5af7f23SKarl Rupp   else {
683bf7f4e0aSPeter Brune     ierr = VecNorm(F, NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
684bf7f4e0aSPeter Brune   }
685bf7f4e0aSPeter Brune 
68657a83faaSBarry Smith   ierr = PetscLogEventBegin(SNESLINESEARCH_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
687bf7f4e0aSPeter Brune 
688bf7f4e0aSPeter Brune   ierr = (*linesearch->ops->apply)(linesearch);CHKERRQ(ierr);
689bf7f4e0aSPeter Brune 
69057a83faaSBarry Smith   ierr = PetscLogEventEnd(SNESLINESEARCH_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
691bf7f4e0aSPeter Brune 
692f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
693bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
694bf7f4e0aSPeter Brune }
695bf7f4e0aSPeter Brune 
696f40b411bSPeter Brune /*@
697f1c6b773SPeter Brune    SNESLineSearchDestroy - Destroys the line search instance.
698f40b411bSPeter Brune 
699f1c6b773SPeter Brune    Collective on SNESLineSearch
700f40b411bSPeter Brune 
701f40b411bSPeter Brune    Input Parameters:
7028e557f58SPeter Brune .  linesearch - The linesearch context
703f40b411bSPeter Brune 
704f40b411bSPeter Brune    Level: Intermediate
705f40b411bSPeter Brune 
70678bcb3b5SPeter Brune .keywords: SNESLineSearch, Destroy
707f40b411bSPeter Brune 
708cd7522eaSPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchReset(), SNESDestroy()
709f40b411bSPeter Brune @*/
710bf388a1fSBarry Smith PetscErrorCode SNESLineSearchDestroy(SNESLineSearch * linesearch)
711bf388a1fSBarry Smith {
712bf7f4e0aSPeter Brune   PetscErrorCode ierr;
713bf388a1fSBarry Smith 
714bf7f4e0aSPeter Brune   PetscFunctionBegin;
715bf7f4e0aSPeter Brune   if (!*linesearch) PetscFunctionReturn(0);
716f1c6b773SPeter Brune   PetscValidHeaderSpecific((*linesearch),SNESLINESEARCH_CLASSID,1);
717bf7f4e0aSPeter Brune   if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = 0; PetscFunctionReturn(0);}
718e04113cfSBarry Smith   ierr = PetscObjectSAWsViewOff((PetscObject)*linesearch);CHKERRQ(ierr);
71922d28d08SBarry Smith   ierr = SNESLineSearchReset(*linesearch);CHKERRQ(ierr);
720f5af7f23SKarl Rupp   if ((*linesearch)->ops->destroy) (*linesearch)->ops->destroy(*linesearch);
721bf7f4e0aSPeter Brune   ierr = PetscViewerDestroy(&(*linesearch)->monitor);CHKERRQ(ierr);
722dcf2fd19SBarry Smith   ierr = SNESLineSearchMonitorCancel((*linesearch));CHKERRQ(ierr);
723e7058c64SPeter Brune   ierr = PetscHeaderDestroy(linesearch);CHKERRQ(ierr);
724bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
725bf7f4e0aSPeter Brune }
726bf7f4e0aSPeter Brune 
727f40b411bSPeter Brune /*@
728dcf2fd19SBarry Smith    SNESLineSearchSetDefaultMonitor - Turns on/off printing useful information and debugging output about the line search.
729bf7f4e0aSPeter Brune 
730bf7f4e0aSPeter Brune    Input Parameters:
731dcf2fd19SBarry Smith +  linesearch - the linesearch object
732dcf2fd19SBarry Smith -  viewer - an ASCII PetscViewer or NULL to turn off monitor
733bf7f4e0aSPeter Brune 
734dcf2fd19SBarry Smith    Logically Collective on SNESLineSearch
735bf7f4e0aSPeter Brune 
736bf7f4e0aSPeter Brune    Options Database:
737dcf2fd19SBarry Smith .   -snes_linesearch_monitor [:filename] - enables the monitor
738bf7f4e0aSPeter Brune 
739bf7f4e0aSPeter Brune    Level: intermediate
740bf7f4e0aSPeter Brune 
741dcf2fd19SBarry Smith    Developer Note: This monitor is implemented differently than the other SNESLineSearchMonitors that are set with
742d12e167eSBarry Smith      SNESLineSearchMonitorSet() since it is called in many locations of the line search routines to display aspects of the
743d12e167eSBarry Smith      line search that are not visible to the other monitors.
744dcf2fd19SBarry Smith 
745d12e167eSBarry Smith .seealso: SNESLineSearchGetDefaultMonitor(), PetscViewer, SNESLineSearchSetMonitor()
746bf7f4e0aSPeter Brune @*/
747dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchSetDefaultMonitor(SNESLineSearch linesearch, PetscViewer viewer)
748bf7f4e0aSPeter Brune {
749bf7f4e0aSPeter Brune   PetscErrorCode ierr;
750bf388a1fSBarry Smith 
751bf7f4e0aSPeter Brune   PetscFunctionBegin;
752dcf2fd19SBarry Smith   if (viewer) {ierr = PetscObjectReference((PetscObject)viewer);CHKERRQ(ierr);}
753bf7f4e0aSPeter Brune   ierr = PetscViewerDestroy(&linesearch->monitor);CHKERRQ(ierr);
754dcf2fd19SBarry Smith   linesearch->monitor = viewer;
755bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
756bf7f4e0aSPeter Brune }
757bf7f4e0aSPeter Brune 
758f40b411bSPeter Brune /*@
759dcf2fd19SBarry Smith    SNESLineSearchGetDefaultMonitor - Gets the PetscViewer instance for the line search monitor.
7606a388c36SPeter Brune 
761f190f2fcSBarry Smith    Input Parameter:
7628e557f58SPeter Brune .  linesearch - linesearch context
763f40b411bSPeter Brune 
764f190f2fcSBarry Smith    Output Parameter:
7658e557f58SPeter Brune .  monitor - monitor context
766f40b411bSPeter Brune 
767f40b411bSPeter Brune    Logically Collective on SNES
768f40b411bSPeter Brune 
769f40b411bSPeter Brune    Options Database Keys:
7708e557f58SPeter Brune .   -snes_linesearch_monitor - enables the monitor
771f40b411bSPeter Brune 
772f40b411bSPeter Brune    Level: intermediate
773f40b411bSPeter Brune 
774dcf2fd19SBarry Smith .seealso: SNESLineSearchSetDefaultMonitor(), PetscViewer
775f40b411bSPeter Brune @*/
776dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchGetDefaultMonitor(SNESLineSearch linesearch, PetscViewer *monitor)
7776a388c36SPeter Brune {
7786a388c36SPeter Brune   PetscFunctionBegin;
779f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
7806a388c36SPeter Brune   if (monitor) {
7816a388c36SPeter Brune     PetscValidPointer(monitor, 2);
7826a388c36SPeter Brune     *monitor = linesearch->monitor;
7836a388c36SPeter Brune   }
7846a388c36SPeter Brune   PetscFunctionReturn(0);
7856a388c36SPeter Brune }
7866a388c36SPeter Brune 
787dcf2fd19SBarry Smith /*@C
788dcf2fd19SBarry Smith    SNESLineSearchMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
789dcf2fd19SBarry Smith 
790dcf2fd19SBarry Smith    Collective on SNESLineSearch
791dcf2fd19SBarry Smith 
792dcf2fd19SBarry Smith    Input Parameters:
793dcf2fd19SBarry Smith +  ls - LineSearch object you wish to monitor
794dcf2fd19SBarry Smith .  name - the monitor type one is seeking
795dcf2fd19SBarry Smith .  help - message indicating what monitoring is done
796dcf2fd19SBarry Smith .  manual - manual page for the monitor
797dcf2fd19SBarry Smith .  monitor - the monitor function
798dcf2fd19SBarry 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
799dcf2fd19SBarry Smith 
800dcf2fd19SBarry Smith    Level: developer
801dcf2fd19SBarry Smith 
802dcf2fd19SBarry Smith .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
803dcf2fd19SBarry Smith           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
804dcf2fd19SBarry Smith           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
805dcf2fd19SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
806dcf2fd19SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
807dcf2fd19SBarry Smith           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
808dcf2fd19SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
809dcf2fd19SBarry Smith @*/
810d12e167eSBarry Smith PetscErrorCode  SNESLineSearchMonitorSetFromOptions(SNESLineSearch ls,const char name[],const char help[], const char manual[],PetscErrorCode (*monitor)(SNESLineSearch,PetscViewerAndFormat*),PetscErrorCode (*monitorsetup)(SNESLineSearch,PetscViewerAndFormat*))
811dcf2fd19SBarry Smith {
812dcf2fd19SBarry Smith   PetscErrorCode    ierr;
813dcf2fd19SBarry Smith   PetscViewer       viewer;
814dcf2fd19SBarry Smith   PetscViewerFormat format;
815dcf2fd19SBarry Smith   PetscBool         flg;
816dcf2fd19SBarry Smith 
817dcf2fd19SBarry Smith   PetscFunctionBegin;
818dcf2fd19SBarry Smith   ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)ls),((PetscObject)ls)->prefix,name,&viewer,&format,&flg);CHKERRQ(ierr);
819dcf2fd19SBarry Smith   if (flg) {
820d12e167eSBarry Smith     PetscViewerAndFormat *vf;
821d12e167eSBarry Smith     ierr = PetscViewerAndFormatCreate(viewer,format,&vf);CHKERRQ(ierr);
822d12e167eSBarry Smith     ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr);
823dcf2fd19SBarry Smith     if (monitorsetup) {
824d12e167eSBarry Smith       ierr = (*monitorsetup)(ls,vf);CHKERRQ(ierr);
825dcf2fd19SBarry Smith     }
826d12e167eSBarry Smith     ierr = SNESLineSearchMonitorSet(ls,(PetscErrorCode (*)(SNESLineSearch,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);CHKERRQ(ierr);
827dcf2fd19SBarry Smith   }
828dcf2fd19SBarry Smith   PetscFunctionReturn(0);
829dcf2fd19SBarry Smith }
830dcf2fd19SBarry Smith 
831f40b411bSPeter Brune /*@
832f1c6b773SPeter Brune    SNESLineSearchSetFromOptions - Sets options for the line search
833f40b411bSPeter Brune 
834f40b411bSPeter Brune    Input Parameters:
8358e557f58SPeter Brune .  linesearch - linesearch context
836f40b411bSPeter Brune 
837f40b411bSPeter Brune    Options Database Keys:
838d4c6564cSPatrick Farrell + -snes_linesearch_type <type> - basic, bt, l2, cp, nleqerr, shell
8393c7d6663SPeter Brune . -snes_linesearch_order <order> - 1, 2, 3.  Most types only support certain orders (bt supports 2 or 3)
8401fe24845SBarry Smith . -snes_linesearch_norms   - Turn on/off the linesearch norms for the basic linesearch typem (SNESLineSearchSetComputeNorms())
84171eef1aeSPeter Brune . -snes_linesearch_minlambda - The minimum step length
8421a9b3a06SPeter Brune . -snes_linesearch_maxstep - The maximum step size
8431a9b3a06SPeter Brune . -snes_linesearch_rtol - Relative tolerance for iterative line searches
8441a9b3a06SPeter Brune . -snes_linesearch_atol - Absolute tolerance for iterative line searches
8451a9b3a06SPeter Brune . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches
8461a9b3a06SPeter Brune . -snes_linesearch_max_it - The number of iterations for iterative line searches
847dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches
848dcf2fd19SBarry Smith . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine
8498e557f58SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter
850cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess.
8511a9b3a06SPeter Brune . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method
8521a9b3a06SPeter Brune - -snes_linesearch_precheck_picard_angle - Angle used in picard precheck method
853f40b411bSPeter Brune 
854f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
855f40b411bSPeter Brune 
856f40b411bSPeter Brune    Level: intermediate
857f40b411bSPeter Brune 
8581fe24845SBarry Smith .seealso: SNESLineSearchCreate(), SNESLineSearchSetOrder(), SNESLineSearchSetType(), SNESLineSearchSetTolerances(), SNESLineSearchSetDamping(), SNESLineSearchPreCheckPicard(),
8591fe24845SBarry Smith           SNESLineSearchType, SNESLineSearchSetComputeNorms()
860f40b411bSPeter Brune @*/
861bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch)
862bf388a1fSBarry Smith {
863bf7f4e0aSPeter Brune   PetscErrorCode    ierr;
8641a4f838cSPeter Brune   const char        *deft = SNESLINESEARCHBASIC;
865bf7f4e0aSPeter Brune   char              type[256];
866bf7f4e0aSPeter Brune   PetscBool         flg, set;
867dcf2fd19SBarry Smith   PetscViewer       viewer;
868bf388a1fSBarry Smith 
869bf7f4e0aSPeter Brune   PetscFunctionBegin;
8700f51fdf8SToby Isaac   ierr = SNESLineSearchRegisterAll();CHKERRQ(ierr);
871bf7f4e0aSPeter Brune 
872bf7f4e0aSPeter Brune   ierr = PetscObjectOptionsBegin((PetscObject)linesearch);CHKERRQ(ierr);
873f5af7f23SKarl Rupp   if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name;
874a264d7a6SBarry Smith   ierr = PetscOptionsFList("-snes_linesearch_type","Linesearch type","SNESLineSearchSetType",SNESLineSearchList,deft,type,256,&flg);CHKERRQ(ierr);
875bf7f4e0aSPeter Brune   if (flg) {
876f1c6b773SPeter Brune     ierr = SNESLineSearchSetType(linesearch,type);CHKERRQ(ierr);
877bf7f4e0aSPeter Brune   } else if (!((PetscObject)linesearch)->type_name) {
878f1c6b773SPeter Brune     ierr = SNESLineSearchSetType(linesearch,deft);CHKERRQ(ierr);
879bf7f4e0aSPeter Brune   }
880bf7f4e0aSPeter Brune 
881dcf2fd19SBarry Smith   ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)linesearch),((PetscObject)linesearch)->prefix,"-snes_linesearch_monitor",&viewer,NULL,&set);CHKERRQ(ierr);
882dcf2fd19SBarry Smith   if (set) {
883dcf2fd19SBarry Smith     ierr = SNESLineSearchSetDefaultMonitor(linesearch,viewer);CHKERRQ(ierr);
884dcf2fd19SBarry Smith     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
885dcf2fd19SBarry Smith   }
886dcf2fd19SBarry Smith   ierr = SNESLineSearchMonitorSetFromOptions(linesearch,"-snes_linesearch_monitor_solution_update","View correction at each iteration","SNESLineSearchMonitorSolutionUpdate",SNESLineSearchMonitorSolutionUpdate,NULL);CHKERRQ(ierr);
887bf7f4e0aSPeter Brune 
8881a9b3a06SPeter Brune   /* tolerances */
88994ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_minlambda","Minimum step length","SNESLineSearchSetTolerances",linesearch->steptol,&linesearch->steptol,NULL);CHKERRQ(ierr);
89094ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_maxstep","Maximum step size","SNESLineSearchSetTolerances",linesearch->maxstep,&linesearch->maxstep,NULL);CHKERRQ(ierr);
89194ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_rtol","Relative tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->rtol,&linesearch->rtol,NULL);CHKERRQ(ierr);
89294ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_atol","Absolute tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->atol,&linesearch->atol,NULL);CHKERRQ(ierr);
89394ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_ltol","Change in lambda tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->ltol,&linesearch->ltol,NULL);CHKERRQ(ierr);
89494ae4db5SBarry Smith   ierr = PetscOptionsInt("-snes_linesearch_max_it","Maximum iterations for iterative line searches","SNESLineSearchSetTolerances",linesearch->max_its,&linesearch->max_its,NULL);CHKERRQ(ierr);
8957a35526eSPeter Brune 
8961a9b3a06SPeter Brune   /* damping parameters */
89794ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_damping","Line search damping and initial step guess","SNESLineSearchSetDamping",linesearch->damping,&linesearch->damping,NULL);CHKERRQ(ierr);
8981a9b3a06SPeter Brune 
89994ae4db5SBarry Smith   ierr = PetscOptionsBool("-snes_linesearch_keeplambda","Use previous lambda as damping","SNESLineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,NULL);CHKERRQ(ierr);
9001a9b3a06SPeter Brune 
9011a9b3a06SPeter Brune   /* precheck */
9027a35526eSPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_precheck_picard","Use a correction that sometimes improves convergence of Picard iteration","SNESLineSearchPreCheckPicard",flg,&flg,&set);CHKERRQ(ierr);
9037a35526eSPeter Brune   if (set) {
9047a35526eSPeter Brune     if (flg) {
9057a35526eSPeter Brune       linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */
906f5af7f23SKarl Rupp 
9077a35526eSPeter Brune       ierr = PetscOptionsReal("-snes_linesearch_precheck_picard_angle","Maximum angle at which to activate the correction",
9080298fd71SBarry Smith                               "none",linesearch->precheck_picard_angle,&linesearch->precheck_picard_angle,NULL);CHKERRQ(ierr);
9097a35526eSPeter Brune       ierr = SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&linesearch->precheck_picard_angle);CHKERRQ(ierr);
9107a35526eSPeter Brune     } else {
9110298fd71SBarry Smith       ierr = SNESLineSearchSetPreCheck(linesearch,NULL,NULL);CHKERRQ(ierr);
9127a35526eSPeter Brune     }
9137a35526eSPeter Brune   }
91494ae4db5SBarry Smith   ierr = PetscOptionsInt("-snes_linesearch_order","Order of approximation used in the line search","SNESLineSearchSetOrder",linesearch->order,&linesearch->order,NULL);CHKERRQ(ierr);
91594ae4db5SBarry Smith   ierr = PetscOptionsBool("-snes_linesearch_norms","Compute final norms in line search","SNESLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,NULL);CHKERRQ(ierr);
9167a35526eSPeter Brune 
9175a9b6599SPeter Brune   if (linesearch->ops->setfromoptions) {
918e55864a3SBarry Smith     (*linesearch->ops->setfromoptions)(PetscOptionsObject,linesearch);CHKERRQ(ierr);
9195a9b6599SPeter Brune   }
9205a9b6599SPeter Brune 
9210633abcbSJed Brown   ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)linesearch);CHKERRQ(ierr);
922bf7f4e0aSPeter Brune   ierr = PetscOptionsEnd();CHKERRQ(ierr);
923bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
924bf7f4e0aSPeter Brune }
925bf7f4e0aSPeter Brune 
926f40b411bSPeter Brune /*@
927f190f2fcSBarry Smith    SNESLineSearchView - Prints useful information about the line search
928f40b411bSPeter Brune 
929f40b411bSPeter Brune    Input Parameters:
9308e557f58SPeter Brune .  linesearch - linesearch context
931f40b411bSPeter Brune 
932f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
933f40b411bSPeter Brune 
934f40b411bSPeter Brune    Level: intermediate
935f40b411bSPeter Brune 
936dcf2fd19SBarry Smith .seealso: SNESLineSearchCreate()
937f40b411bSPeter Brune @*/
938bf388a1fSBarry Smith PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer)
939bf388a1fSBarry Smith {
9407f1410a3SPeter Brune   PetscErrorCode ierr;
9417f1410a3SPeter Brune   PetscBool      iascii;
942bf388a1fSBarry Smith 
943bf7f4e0aSPeter Brune   PetscFunctionBegin;
9447f1410a3SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
9457f1410a3SPeter Brune   if (!viewer) {
946ce94432eSBarry Smith     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch),&viewer);CHKERRQ(ierr);
9477f1410a3SPeter Brune   }
9487f1410a3SPeter Brune   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
9497f1410a3SPeter Brune   PetscCheckSameComm(linesearch,1,viewer,2);
950f40b411bSPeter Brune 
951251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
9527f1410a3SPeter Brune   if (iascii) {
953dae58748SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)linesearch,viewer);CHKERRQ(ierr);
9547f1410a3SPeter Brune     if (linesearch->ops->view) {
9557f1410a3SPeter Brune       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
9567f1410a3SPeter Brune       ierr = (*linesearch->ops->view)(linesearch,viewer);CHKERRQ(ierr);
9577f1410a3SPeter Brune       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
9587f1410a3SPeter Brune     }
959c69d1a72SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep,(double)linesearch->steptol);CHKERRQ(ierr);
960c69d1a72SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol,(double)linesearch->atol,(double)linesearch->ltol);CHKERRQ(ierr);
9617f1410a3SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D\n", linesearch->max_its);CHKERRQ(ierr);
9626b2b7091SBarry Smith     if (linesearch->ops->precheck) {
9636b2b7091SBarry Smith       if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) {
9647f1410a3SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  using precheck step to speed up Picard convergence\n", linesearch->max_its);CHKERRQ(ierr);
9657f1410a3SPeter Brune       } else {
9667f1410a3SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  using user-defined precheck step\n", linesearch->max_its);CHKERRQ(ierr);
9677f1410a3SPeter Brune       }
9687f1410a3SPeter Brune     }
9696b2b7091SBarry Smith     if (linesearch->ops->postcheck) {
9707f1410a3SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer,"  using user-defined postcheck step\n", linesearch->max_its);CHKERRQ(ierr);
9717f1410a3SPeter Brune     }
9727f1410a3SPeter Brune   }
973bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
974bf7f4e0aSPeter Brune }
975bf7f4e0aSPeter Brune 
976ea5d4fccSPeter Brune /*@C
977f1c6b773SPeter Brune    SNESLineSearchSetType - Sets the linesearch type
978f40b411bSPeter Brune 
979f190f2fcSBarry Smith    Logically Collective on SNESLineSearch
980f190f2fcSBarry Smith 
981f40b411bSPeter Brune    Input Parameters:
9828e557f58SPeter Brune +  linesearch - linesearch context
983f40b411bSPeter Brune -  type - The type of line search to be used
984f40b411bSPeter Brune 
985cd7522eaSPeter Brune    Available Types:
9861fe24845SBarry Smith +  SNESLINESEARCHBASIC - Simple damping line search, defaults to using the full Newton step
9871fe24845SBarry Smith .  SNESLINESEARCHBT - Backtracking line search over the L2 norm of the function
9881fe24845SBarry Smith .  SNESLINESEARCHL2 - Secant line search over the L2 norm of the function
9891fe24845SBarry Smith .  SNESLINESEARCHCP - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x)
9901fe24845SBarry Smith .  SNESLINESEARCHNLEQERR - Affine-covariant error-oriented linesearch
9911fe24845SBarry Smith -  SNESLINESEARCHSHELL - User provided SNESLineSearch implementation
9921fe24845SBarry Smith 
9931fe24845SBarry Smith    Options Database:
9941fe24845SBarry Smith .  -snes_linesearch_type <type> - basic, bt, l2, cp, nleqerr, shell
995cd7522eaSPeter Brune 
996f40b411bSPeter Brune    Level: intermediate
997f40b411bSPeter Brune 
9981fe24845SBarry Smith .seealso: SNESLineSearchCreate(), SNESLineSearchType, SNESLineSearchSetFromOptions()
999f40b411bSPeter Brune @*/
100019fd82e9SBarry Smith PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type)
1001bf7f4e0aSPeter Brune {
1002f1c6b773SPeter Brune   PetscErrorCode ierr,(*r)(SNESLineSearch);
1003bf7f4e0aSPeter Brune   PetscBool      match;
1004bf7f4e0aSPeter Brune 
1005bf7f4e0aSPeter Brune   PetscFunctionBegin;
1006f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1007bf7f4e0aSPeter Brune   PetscValidCharPointer(type,2);
1008bf7f4e0aSPeter Brune 
1009251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)linesearch,type,&match);CHKERRQ(ierr);
1010bf7f4e0aSPeter Brune   if (match) PetscFunctionReturn(0);
1011bf7f4e0aSPeter Brune 
10121c9cd337SJed Brown   ierr = PetscFunctionListFind(SNESLineSearchList,type,&r);CHKERRQ(ierr);
1013bf7f4e0aSPeter Brune   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type);
1014bf7f4e0aSPeter Brune   /* Destroy the previous private linesearch context */
1015bf7f4e0aSPeter Brune   if (linesearch->ops->destroy) {
1016bf7f4e0aSPeter Brune     ierr = (*(linesearch)->ops->destroy)(linesearch);CHKERRQ(ierr);
1017f5af7f23SKarl Rupp 
10180298fd71SBarry Smith     linesearch->ops->destroy = NULL;
1019bf7f4e0aSPeter Brune   }
1020f1c6b773SPeter Brune   /* Reinitialize function pointers in SNESLineSearchOps structure */
1021bf7f4e0aSPeter Brune   linesearch->ops->apply          = 0;
1022bf7f4e0aSPeter Brune   linesearch->ops->view           = 0;
1023bf7f4e0aSPeter Brune   linesearch->ops->setfromoptions = 0;
1024bf7f4e0aSPeter Brune   linesearch->ops->destroy        = 0;
1025bf7f4e0aSPeter Brune 
1026bf7f4e0aSPeter Brune   ierr = PetscObjectChangeTypeName((PetscObject)linesearch,type);CHKERRQ(ierr);
1027bf7f4e0aSPeter Brune   ierr = (*r)(linesearch);CHKERRQ(ierr);
1028bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1029bf7f4e0aSPeter Brune }
1030bf7f4e0aSPeter Brune 
1031f40b411bSPeter Brune /*@
103278bcb3b5SPeter Brune    SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation.
1033f40b411bSPeter Brune 
1034f40b411bSPeter Brune    Input Parameters:
10358e557f58SPeter Brune +  linesearch - linesearch context
1036f40b411bSPeter Brune -  snes - The snes instance
1037f40b411bSPeter Brune 
103878bcb3b5SPeter Brune    Level: developer
103978bcb3b5SPeter Brune 
104078bcb3b5SPeter Brune    Notes:
1041f190f2fcSBarry Smith    This happens automatically when the line search is obtained/created with
10427601faf0SJed Brown    SNESGetLineSearch().  This routine is therefore mainly called within SNES
104378bcb3b5SPeter Brune    implementations.
1044f40b411bSPeter Brune 
10458141a3b9SPeter Brune    Level: developer
1046f40b411bSPeter Brune 
1047cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
1048f40b411bSPeter Brune @*/
1049bf388a1fSBarry Smith PetscErrorCode  SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes)
1050bf388a1fSBarry Smith {
1051bf7f4e0aSPeter Brune   PetscFunctionBegin;
1052f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1053bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,2);
1054bf7f4e0aSPeter Brune   linesearch->snes = snes;
1055bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1056bf7f4e0aSPeter Brune }
1057bf7f4e0aSPeter Brune 
1058f40b411bSPeter Brune /*@
10598141a3b9SPeter Brune    SNESLineSearchGetSNES - Gets the SNES instance associated with the line search.
10608141a3b9SPeter Brune    Having an associated SNES is necessary because most line search implementations must be able to
10618141a3b9SPeter Brune    evaluate the function using SNESComputeFunction() for the associated SNES.  This routine
10628141a3b9SPeter Brune    is used in the line search implementations when one must get this associated SNES instance.
1063f40b411bSPeter Brune 
1064f40b411bSPeter Brune    Input Parameters:
10658e557f58SPeter Brune .  linesearch - linesearch context
1066f40b411bSPeter Brune 
1067f40b411bSPeter Brune    Output Parameters:
1068f40b411bSPeter Brune .  snes - The snes instance
1069f40b411bSPeter Brune 
10708141a3b9SPeter Brune    Level: developer
1071f40b411bSPeter Brune 
1072cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
1073f40b411bSPeter Brune @*/
1074bf388a1fSBarry Smith PetscErrorCode  SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes)
1075bf388a1fSBarry Smith {
1076bf7f4e0aSPeter Brune   PetscFunctionBegin;
1077f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
10786a388c36SPeter Brune   PetscValidPointer(snes, 2);
1079bf7f4e0aSPeter Brune   *snes = linesearch->snes;
1080bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1081bf7f4e0aSPeter Brune }
1082bf7f4e0aSPeter Brune 
1083f40b411bSPeter Brune /*@
1084f1c6b773SPeter Brune    SNESLineSearchGetLambda - Gets the last linesearch steplength discovered.
1085f40b411bSPeter Brune 
1086f40b411bSPeter Brune    Input Parameters:
10878e557f58SPeter Brune .  linesearch - linesearch context
1088f40b411bSPeter Brune 
1089f40b411bSPeter Brune    Output Parameters:
1090cd7522eaSPeter Brune .  lambda - The last steplength computed during SNESLineSearchApply()
1091f40b411bSPeter Brune 
109278bcb3b5SPeter Brune    Level: advanced
109378bcb3b5SPeter Brune 
10948e557f58SPeter Brune    Notes:
10958e557f58SPeter Brune    This is useful in methods where the solver is ill-scaled and
109678bcb3b5SPeter Brune    requires some adaptive notion of the difference in scale between the
109778bcb3b5SPeter Brune    solution and the function.  For instance, SNESQN may be scaled by the
109878bcb3b5SPeter Brune    line search lambda using the argument -snes_qn_scaling ls.
109978bcb3b5SPeter Brune 
1100cd7522eaSPeter Brune .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetDamping(), SNESLineSearchApply()
1101f40b411bSPeter Brune @*/
1102f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda)
11036a388c36SPeter Brune {
11046a388c36SPeter Brune   PetscFunctionBegin;
1105f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
11066a388c36SPeter Brune   PetscValidPointer(lambda, 2);
11076a388c36SPeter Brune   *lambda = linesearch->lambda;
11086a388c36SPeter Brune   PetscFunctionReturn(0);
11096a388c36SPeter Brune }
11106a388c36SPeter Brune 
1111f40b411bSPeter Brune /*@
1112f1c6b773SPeter Brune    SNESLineSearchSetLambda - Sets the linesearch steplength.
1113f40b411bSPeter Brune 
1114f40b411bSPeter Brune    Input Parameters:
11158e557f58SPeter Brune +  linesearch - linesearch context
1116f40b411bSPeter Brune -  lambda - The last steplength.
1117f40b411bSPeter Brune 
1118cd7522eaSPeter Brune    Notes:
1119f190f2fcSBarry Smith    This routine is typically used within implementations of SNESLineSearchApply()
1120cd7522eaSPeter Brune    to set the final steplength.  This routine (and SNESLineSearchGetLambda()) were
1121cd7522eaSPeter Brune    added in order to facilitate Quasi-Newton methods that use the previous steplength
1122cd7522eaSPeter Brune    as an inner scaling parameter.
1123cd7522eaSPeter Brune 
112478bcb3b5SPeter Brune    Level: advanced
1125f40b411bSPeter Brune 
1126f1c6b773SPeter Brune .seealso: SNESLineSearchGetLambda()
1127f40b411bSPeter Brune @*/
1128f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
11296a388c36SPeter Brune {
11306a388c36SPeter Brune   PetscFunctionBegin;
1131f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
11326a388c36SPeter Brune   linesearch->lambda = lambda;
11336a388c36SPeter Brune   PetscFunctionReturn(0);
11346a388c36SPeter Brune }
11356a388c36SPeter Brune 
1136f40b411bSPeter Brune /*@
11373c7d6663SPeter Brune    SNESLineSearchGetTolerances - Gets the tolerances for the linesearch.  These include
113878bcb3b5SPeter Brune    tolerances for the relative and absolute change in the function norm, the change
113978bcb3b5SPeter Brune    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
114078bcb3b5SPeter Brune    and the maximum number of iterations the line search procedure may take.
1141f40b411bSPeter Brune 
1142f40b411bSPeter Brune    Input Parameters:
11438e557f58SPeter Brune .  linesearch - linesearch context
1144f40b411bSPeter Brune 
1145f40b411bSPeter Brune    Output Parameters:
1146516fe3c3SPeter Brune +  steptol - The minimum steplength
11476cc8e53bSPeter Brune .  maxstep - The maximum steplength
1148516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
1149516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
1150516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
1151516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
1152f40b411bSPeter Brune 
115378bcb3b5SPeter Brune    Level: intermediate
115478bcb3b5SPeter Brune 
115578bcb3b5SPeter Brune    Notes:
115678bcb3b5SPeter Brune    Different line searches may implement these parameters slightly differently as
11573c7d6663SPeter Brune    the type requires.
1158516fe3c3SPeter Brune 
1159f1c6b773SPeter Brune .seealso: SNESLineSearchSetTolerances()
1160f40b411bSPeter Brune @*/
1161f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
11626a388c36SPeter Brune {
11636a388c36SPeter Brune   PetscFunctionBegin;
1164f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1165516fe3c3SPeter Brune   if (steptol) {
11666a388c36SPeter Brune     PetscValidPointer(steptol, 2);
11676a388c36SPeter Brune     *steptol = linesearch->steptol;
1168516fe3c3SPeter Brune   }
1169516fe3c3SPeter Brune   if (maxstep) {
1170516fe3c3SPeter Brune     PetscValidPointer(maxstep, 3);
1171516fe3c3SPeter Brune     *maxstep = linesearch->maxstep;
1172516fe3c3SPeter Brune   }
1173516fe3c3SPeter Brune   if (rtol) {
1174516fe3c3SPeter Brune     PetscValidPointer(rtol, 4);
1175516fe3c3SPeter Brune     *rtol = linesearch->rtol;
1176516fe3c3SPeter Brune   }
1177516fe3c3SPeter Brune   if (atol) {
1178516fe3c3SPeter Brune     PetscValidPointer(atol, 5);
1179516fe3c3SPeter Brune     *atol = linesearch->atol;
1180516fe3c3SPeter Brune   }
1181516fe3c3SPeter Brune   if (ltol) {
1182516fe3c3SPeter Brune     PetscValidPointer(ltol, 6);
1183516fe3c3SPeter Brune     *ltol = linesearch->ltol;
1184516fe3c3SPeter Brune   }
1185516fe3c3SPeter Brune   if (max_its) {
1186516fe3c3SPeter Brune     PetscValidPointer(max_its, 7);
1187516fe3c3SPeter Brune     *max_its = linesearch->max_its;
1188516fe3c3SPeter Brune   }
11896a388c36SPeter Brune   PetscFunctionReturn(0);
11906a388c36SPeter Brune }
11916a388c36SPeter Brune 
1192f40b411bSPeter Brune /*@
11933c7d6663SPeter Brune    SNESLineSearchSetTolerances -  Gets the tolerances for the linesearch.  These include
119478bcb3b5SPeter Brune    tolerances for the relative and absolute change in the function norm, the change
119578bcb3b5SPeter Brune    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
119678bcb3b5SPeter Brune    and the maximum number of iterations the line search procedure may take.
1197f40b411bSPeter Brune 
1198f40b411bSPeter Brune    Input Parameters:
11998e557f58SPeter Brune +  linesearch - linesearch context
1200516fe3c3SPeter Brune .  steptol - The minimum steplength
12016cc8e53bSPeter Brune .  maxstep - The maximum steplength
1202516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
1203516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
1204516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
1205516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
1206f40b411bSPeter Brune 
120778bcb3b5SPeter Brune    Notes:
12083c7d6663SPeter Brune    The user may choose to not set any of the tolerances using PETSC_DEFAULT in
120978bcb3b5SPeter Brune    place of an argument.
1210f40b411bSPeter Brune 
121178bcb3b5SPeter Brune    Level: intermediate
1212516fe3c3SPeter Brune 
1213f1c6b773SPeter Brune .seealso: SNESLineSearchGetTolerances()
1214f40b411bSPeter Brune @*/
1215f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its)
12166a388c36SPeter Brune {
12176a388c36SPeter Brune   PetscFunctionBegin;
1218f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1219d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,steptol,2);
1220d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,maxstep,3);
1221d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,rtol,4);
1222d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,atol,5);
1223d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,ltol,6);
1224d3952184SSatish Balay   PetscValidLogicalCollectiveInt(linesearch,max_its,7);
1225d3952184SSatish Balay 
1226d3952184SSatish Balay   if (steptol!= PETSC_DEFAULT) {
1227ce94432eSBarry Smith     if (steptol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Minimum step length %14.12e must be non-negative",(double)steptol);
12286a388c36SPeter Brune     linesearch->steptol = steptol;
1229d3952184SSatish Balay   }
1230d3952184SSatish Balay 
1231d3952184SSatish Balay   if (maxstep!= PETSC_DEFAULT) {
1232ce94432eSBarry Smith     if (maxstep < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum step length %14.12e must be non-negative",(double)maxstep);
1233516fe3c3SPeter Brune     linesearch->maxstep = maxstep;
1234d3952184SSatish Balay   }
1235d3952184SSatish Balay 
1236d3952184SSatish Balay   if (rtol != PETSC_DEFAULT) {
1237ce94432eSBarry 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);
1238516fe3c3SPeter Brune     linesearch->rtol = rtol;
1239d3952184SSatish Balay   }
1240d3952184SSatish Balay 
1241d3952184SSatish Balay   if (atol != PETSC_DEFAULT) {
1242ce94432eSBarry Smith     if (atol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %14.12e must be non-negative",(double)atol);
1243516fe3c3SPeter Brune     linesearch->atol = atol;
1244d3952184SSatish Balay   }
1245d3952184SSatish Balay 
1246d3952184SSatish Balay   if (ltol != PETSC_DEFAULT) {
1247ce94432eSBarry Smith     if (ltol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Labmda tolerance %14.12e must be non-negative",(double)ltol);
1248516fe3c3SPeter Brune     linesearch->ltol = ltol;
1249d3952184SSatish Balay   }
1250d3952184SSatish Balay 
1251d3952184SSatish Balay   if (max_its != PETSC_DEFAULT) {
1252ce94432eSBarry Smith     if (max_its < 0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",max_its);
1253516fe3c3SPeter Brune     linesearch->max_its = max_its;
1254d3952184SSatish Balay   }
12556a388c36SPeter Brune   PetscFunctionReturn(0);
12566a388c36SPeter Brune }
12576a388c36SPeter Brune 
1258f40b411bSPeter Brune /*@
1259f1c6b773SPeter Brune    SNESLineSearchGetDamping - Gets the line search damping parameter.
1260f40b411bSPeter Brune 
1261f40b411bSPeter Brune    Input Parameters:
12628e557f58SPeter Brune .  linesearch - linesearch context
1263f40b411bSPeter Brune 
1264f40b411bSPeter Brune    Output Parameters:
12658e557f58SPeter Brune .  damping - The damping parameter
1266f40b411bSPeter Brune 
126778bcb3b5SPeter Brune    Level: advanced
1268f40b411bSPeter Brune 
126978bcb3b5SPeter Brune .seealso: SNESLineSearchGetStepTolerance(), SNESQN
1270f40b411bSPeter Brune @*/
1271f40b411bSPeter Brune 
1272f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping)
12736a388c36SPeter Brune {
12746a388c36SPeter Brune   PetscFunctionBegin;
1275f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
12766a388c36SPeter Brune   PetscValidPointer(damping, 2);
12776a388c36SPeter Brune   *damping = linesearch->damping;
12786a388c36SPeter Brune   PetscFunctionReturn(0);
12796a388c36SPeter Brune }
12806a388c36SPeter Brune 
1281f40b411bSPeter Brune /*@
1282f1c6b773SPeter Brune    SNESLineSearchSetDamping - Sets the line search damping paramter.
1283f40b411bSPeter Brune 
1284f40b411bSPeter Brune    Input Parameters:
128503fd4120SBarry Smith +  linesearch - linesearch context
128603fd4120SBarry Smith -  damping - The damping parameter
1287f40b411bSPeter Brune 
128803fd4120SBarry Smith    Options Database:
128903fd4120SBarry Smith .   -snes_linesearch_damping
1290f40b411bSPeter Brune    Level: intermediate
1291f40b411bSPeter Brune 
1292cd7522eaSPeter Brune    Notes:
1293cd7522eaSPeter Brune    The basic line search merely takes the update step scaled by the damping parameter.
1294cd7522eaSPeter Brune    The use of the damping parameter in the l2 and cp line searches is much more subtle;
129578bcb3b5SPeter Brune    it is used as a starting point in calculating the secant step. However, the eventual
1296cd7522eaSPeter Brune    step may be of greater length than the damping parameter.  In the bt line search it is
1297cd7522eaSPeter Brune    used as the maximum possible step length, as the bt line search only backtracks.
1298cd7522eaSPeter Brune 
1299f1c6b773SPeter Brune .seealso: SNESLineSearchGetDamping()
1300f40b411bSPeter Brune @*/
1301f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping)
13026a388c36SPeter Brune {
13036a388c36SPeter Brune   PetscFunctionBegin;
1304f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
13056a388c36SPeter Brune   linesearch->damping = damping;
13066a388c36SPeter Brune   PetscFunctionReturn(0);
13076a388c36SPeter Brune }
13086a388c36SPeter Brune 
130959405d5eSPeter Brune /*@
131059405d5eSPeter Brune    SNESLineSearchGetOrder - Gets the line search approximation order.
131159405d5eSPeter Brune 
131259405d5eSPeter Brune    Input Parameters:
131378bcb3b5SPeter Brune .  linesearch - linesearch context
131459405d5eSPeter Brune 
131559405d5eSPeter Brune    Output Parameters:
13168e557f58SPeter Brune .  order - The order
131759405d5eSPeter Brune 
131878bcb3b5SPeter Brune    Possible Values for order:
13193c7d6663SPeter Brune +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
13203c7d6663SPeter Brune .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
13213c7d6663SPeter Brune -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
132278bcb3b5SPeter Brune 
132359405d5eSPeter Brune    Level: intermediate
132459405d5eSPeter Brune 
132559405d5eSPeter Brune .seealso: SNESLineSearchSetOrder()
132659405d5eSPeter Brune @*/
132759405d5eSPeter Brune 
1328b000cd8dSPeter Brune PetscErrorCode  SNESLineSearchGetOrder(SNESLineSearch linesearch,PetscInt *order)
132959405d5eSPeter Brune {
133059405d5eSPeter Brune   PetscFunctionBegin;
133159405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
133259405d5eSPeter Brune   PetscValidPointer(order, 2);
133359405d5eSPeter Brune   *order = linesearch->order;
133459405d5eSPeter Brune   PetscFunctionReturn(0);
133559405d5eSPeter Brune }
133659405d5eSPeter Brune 
133759405d5eSPeter Brune /*@
1338*1f8196a2SJed Brown    SNESLineSearchSetOrder - Sets the maximum order of the polynomial fit used in the line search
133959405d5eSPeter Brune 
134059405d5eSPeter Brune    Input Parameters:
134178bcb3b5SPeter Brune .  linesearch - linesearch context
134278bcb3b5SPeter Brune .  order - The damping parameter
134359405d5eSPeter Brune 
134459405d5eSPeter Brune    Level: intermediate
134559405d5eSPeter Brune 
134678bcb3b5SPeter Brune    Possible Values for order:
13473c7d6663SPeter Brune +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
13483c7d6663SPeter Brune .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
13493c7d6663SPeter Brune -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
135078bcb3b5SPeter Brune 
135159405d5eSPeter Brune    Notes:
135259405d5eSPeter Brune    Variable orders are supported by the following line searches:
135378bcb3b5SPeter Brune +  bt - cubic and quadratic
135478bcb3b5SPeter Brune -  cp - linear and quadratic
135559405d5eSPeter Brune 
1356*1f8196a2SJed Brown .seealso: SNESLineSearchGetOrder(), SNESLineSearchSetDamping()
135759405d5eSPeter Brune @*/
1358b000cd8dSPeter Brune PetscErrorCode  SNESLineSearchSetOrder(SNESLineSearch linesearch,PetscInt order)
135959405d5eSPeter Brune {
136059405d5eSPeter Brune   PetscFunctionBegin;
136159405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
136259405d5eSPeter Brune   linesearch->order = order;
136359405d5eSPeter Brune   PetscFunctionReturn(0);
136459405d5eSPeter Brune }
136559405d5eSPeter Brune 
1366f40b411bSPeter Brune /*@
1367f1c6b773SPeter Brune    SNESLineSearchGetNorms - Gets the norms for for X, Y, and F.
1368f40b411bSPeter Brune 
1369f40b411bSPeter Brune    Input Parameters:
137078bcb3b5SPeter Brune .  linesearch - linesearch context
1371f40b411bSPeter Brune 
1372f40b411bSPeter Brune    Output Parameters:
1373f40b411bSPeter Brune +  xnorm - The norm of the current solution
1374f40b411bSPeter Brune .  fnorm - The norm of the current function
1375f40b411bSPeter Brune -  ynorm - The norm of the current update
1376f40b411bSPeter Brune 
1377cd7522eaSPeter Brune    Notes:
1378cd7522eaSPeter Brune    This function is mainly called from SNES implementations.
1379cd7522eaSPeter Brune 
138078bcb3b5SPeter Brune    Level: developer
1381f40b411bSPeter Brune 
1382f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms() SNESLineSearchGetVecs()
1383f40b411bSPeter Brune @*/
1384f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm)
1385bf7f4e0aSPeter Brune {
1386bf7f4e0aSPeter Brune   PetscFunctionBegin;
1387f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1388f5af7f23SKarl Rupp   if (xnorm) *xnorm = linesearch->xnorm;
1389f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
1390f5af7f23SKarl Rupp   if (ynorm) *ynorm = linesearch->ynorm;
1391bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1392bf7f4e0aSPeter Brune }
1393bf7f4e0aSPeter Brune 
1394f40b411bSPeter Brune /*@
1395f1c6b773SPeter Brune    SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F.
1396f40b411bSPeter Brune 
1397f40b411bSPeter Brune    Input Parameters:
139878bcb3b5SPeter Brune +  linesearch - linesearch context
1399f40b411bSPeter Brune .  xnorm - The norm of the current solution
1400f40b411bSPeter Brune .  fnorm - The norm of the current function
1401f40b411bSPeter Brune -  ynorm - The norm of the current update
1402f40b411bSPeter Brune 
140378bcb3b5SPeter Brune    Level: advanced
1404f40b411bSPeter Brune 
1405f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1406f40b411bSPeter Brune @*/
1407f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
14086a388c36SPeter Brune {
14096a388c36SPeter Brune   PetscFunctionBegin;
1410f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
14116a388c36SPeter Brune   linesearch->xnorm = xnorm;
14126a388c36SPeter Brune   linesearch->fnorm = fnorm;
14136a388c36SPeter Brune   linesearch->ynorm = ynorm;
14146a388c36SPeter Brune   PetscFunctionReturn(0);
14156a388c36SPeter Brune }
14166a388c36SPeter Brune 
1417f40b411bSPeter Brune /*@
1418f1c6b773SPeter Brune    SNESLineSearchComputeNorms - Computes the norms of X, F, and Y.
1419f40b411bSPeter Brune 
1420f40b411bSPeter Brune    Input Parameters:
142178bcb3b5SPeter Brune .  linesearch - linesearch context
1422f40b411bSPeter Brune 
1423f40b411bSPeter Brune    Options Database Keys:
14248e557f58SPeter Brune .   -snes_linesearch_norms - turn norm computation on or off
1425f40b411bSPeter Brune 
1426f40b411bSPeter Brune    Level: intermediate
1427f40b411bSPeter Brune 
142878bcb3b5SPeter Brune .seealso: SNESLineSearchGetNorms, SNESLineSearchSetNorms(), SNESLineSearchSetComputeNorms()
1429f40b411bSPeter Brune @*/
1430f1c6b773SPeter Brune PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
14316a388c36SPeter Brune {
14326a388c36SPeter Brune   PetscErrorCode ierr;
14339bd66eb0SPeter Brune   SNES           snes;
1434bf388a1fSBarry Smith 
14356a388c36SPeter Brune   PetscFunctionBegin;
14366a388c36SPeter Brune   if (linesearch->norms) {
14379bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
1438f1c6b773SPeter Brune       ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
14399bd66eb0SPeter Brune       ierr = VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
14409bd66eb0SPeter Brune       ierr = VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
14419bd66eb0SPeter Brune       ierr = (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);CHKERRQ(ierr);
14429bd66eb0SPeter Brune     } else {
14436a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_func,   NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
14446a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_sol,    NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
14456a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
14466a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_func,     NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
14476a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_sol,      NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
14486a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_update,   NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
14496a388c36SPeter Brune     }
14509bd66eb0SPeter Brune   }
14516a388c36SPeter Brune   PetscFunctionReturn(0);
14526a388c36SPeter Brune }
14536a388c36SPeter Brune 
14546f263ca3SPeter Brune /*@
14556f263ca3SPeter Brune    SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search.
14566f263ca3SPeter Brune 
14576f263ca3SPeter Brune    Input Parameters:
145878bcb3b5SPeter Brune +  linesearch  - linesearch context
145978bcb3b5SPeter Brune -  flg  - indicates whether or not to compute norms
14606f263ca3SPeter Brune 
14616f263ca3SPeter Brune    Options Database Keys:
1462*1f8196a2SJed Brown .   -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic linesearch
14636f263ca3SPeter Brune 
14646f263ca3SPeter Brune    Notes:
1465*1f8196a2SJed Brown    This is most relevant to the SNESLINESEARCHBASIC line search type since most line searches have a stopping criteria involving the norm.
14666f263ca3SPeter Brune 
14676f263ca3SPeter Brune    Level: intermediate
14686f263ca3SPeter Brune 
14691a4f838cSPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetNorms(), SNESLineSearchComputeNorms(), SNESLINESEARCHBASIC
14706f263ca3SPeter Brune @*/
14716f263ca3SPeter Brune PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg)
14726f263ca3SPeter Brune {
14736f263ca3SPeter Brune   PetscFunctionBegin;
14746f263ca3SPeter Brune   linesearch->norms = flg;
14756f263ca3SPeter Brune   PetscFunctionReturn(0);
14766f263ca3SPeter Brune }
14776f263ca3SPeter Brune 
1478f40b411bSPeter Brune /*@
1479f1c6b773SPeter Brune    SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context
1480f40b411bSPeter Brune 
1481f40b411bSPeter Brune    Input Parameters:
148278bcb3b5SPeter Brune .  linesearch - linesearch context
1483f40b411bSPeter Brune 
1484f40b411bSPeter Brune    Output Parameters:
14856232e825SPeter Brune +  X - Solution vector
14866232e825SPeter Brune .  F - Function vector
14876232e825SPeter Brune .  Y - Search direction vector
14886232e825SPeter Brune .  W - Solution work vector
14896232e825SPeter Brune -  G - Function work vector
14906232e825SPeter Brune 
14917bba9028SPeter Brune    Notes:
14927bba9028SPeter Brune    At the beginning of a line search application, X should contain a
14936232e825SPeter Brune    solution and the vector F the function computed at X.  At the end of the
14946232e825SPeter Brune    line search application, X should contain the new solution, and F the
14956232e825SPeter Brune    function evaluated at the new solution.
1496f40b411bSPeter Brune 
14972a7a6963SBarry Smith    These vectors are owned by the SNESLineSearch and should not be destroyed by the caller
14982a7a6963SBarry Smith 
149978bcb3b5SPeter Brune    Level: advanced
1500f40b411bSPeter Brune 
1501f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1502f40b411bSPeter Brune @*/
1503bf388a1fSBarry Smith PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G)
1504bf388a1fSBarry Smith {
15056a388c36SPeter Brune   PetscFunctionBegin;
1506f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
15076a388c36SPeter Brune   if (X) {
15086a388c36SPeter Brune     PetscValidPointer(X, 2);
15096a388c36SPeter Brune     *X = linesearch->vec_sol;
15106a388c36SPeter Brune   }
15116a388c36SPeter Brune   if (F) {
15126a388c36SPeter Brune     PetscValidPointer(F, 3);
15136a388c36SPeter Brune     *F = linesearch->vec_func;
15146a388c36SPeter Brune   }
15156a388c36SPeter Brune   if (Y) {
15166a388c36SPeter Brune     PetscValidPointer(Y, 4);
15176a388c36SPeter Brune     *Y = linesearch->vec_update;
15186a388c36SPeter Brune   }
15196a388c36SPeter Brune   if (W) {
15206a388c36SPeter Brune     PetscValidPointer(W, 5);
15216a388c36SPeter Brune     *W = linesearch->vec_sol_new;
15226a388c36SPeter Brune   }
15236a388c36SPeter Brune   if (G) {
15246a388c36SPeter Brune     PetscValidPointer(G, 6);
15256a388c36SPeter Brune     *G = linesearch->vec_func_new;
15266a388c36SPeter Brune   }
15276a388c36SPeter Brune   PetscFunctionReturn(0);
15286a388c36SPeter Brune }
15296a388c36SPeter Brune 
1530f40b411bSPeter Brune /*@
1531f1c6b773SPeter Brune    SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context
1532f40b411bSPeter Brune 
1533f40b411bSPeter Brune    Input Parameters:
153478bcb3b5SPeter Brune +  linesearch - linesearch context
15356232e825SPeter Brune .  X - Solution vector
15366232e825SPeter Brune .  F - Function vector
15376232e825SPeter Brune .  Y - Search direction vector
15386232e825SPeter Brune .  W - Solution work vector
15396232e825SPeter Brune -  G - Function work vector
1540f40b411bSPeter Brune 
154178bcb3b5SPeter Brune    Level: advanced
1542f40b411bSPeter Brune 
1543f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms(), SNESLineSearchGetVecs()
1544f40b411bSPeter Brune @*/
1545bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G)
1546bf388a1fSBarry Smith {
15476a388c36SPeter Brune   PetscFunctionBegin;
1548f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
15496a388c36SPeter Brune   if (X) {
15506a388c36SPeter Brune     PetscValidHeaderSpecific(X,VEC_CLASSID,2);
15516a388c36SPeter Brune     linesearch->vec_sol = X;
15526a388c36SPeter Brune   }
15536a388c36SPeter Brune   if (F) {
15546a388c36SPeter Brune     PetscValidHeaderSpecific(F,VEC_CLASSID,3);
15556a388c36SPeter Brune     linesearch->vec_func = F;
15566a388c36SPeter Brune   }
15576a388c36SPeter Brune   if (Y) {
15586a388c36SPeter Brune     PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
15596a388c36SPeter Brune     linesearch->vec_update = Y;
15606a388c36SPeter Brune   }
15616a388c36SPeter Brune   if (W) {
15626a388c36SPeter Brune     PetscValidHeaderSpecific(W,VEC_CLASSID,5);
15636a388c36SPeter Brune     linesearch->vec_sol_new = W;
15646a388c36SPeter Brune   }
15656a388c36SPeter Brune   if (G) {
15666a388c36SPeter Brune     PetscValidHeaderSpecific(G,VEC_CLASSID,6);
15676a388c36SPeter Brune     linesearch->vec_func_new = G;
15686a388c36SPeter Brune   }
15696a388c36SPeter Brune   PetscFunctionReturn(0);
15706a388c36SPeter Brune }
15716a388c36SPeter Brune 
1572e7058c64SPeter Brune /*@C
1573f1c6b773SPeter Brune    SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1574e7058c64SPeter Brune    SNES options in the database.
1575e7058c64SPeter Brune 
1576cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
1577e7058c64SPeter Brune 
1578e7058c64SPeter Brune    Input Parameters:
1579e7058c64SPeter Brune +  snes - the SNES context
1580e7058c64SPeter Brune -  prefix - the prefix to prepend to all option names
1581e7058c64SPeter Brune 
1582e7058c64SPeter Brune    Notes:
1583e7058c64SPeter Brune    A hyphen (-) must NOT be given at the beginning of the prefix name.
1584e7058c64SPeter Brune    The first character of all runtime options is AUTOMATICALLY the hyphen.
1585e7058c64SPeter Brune 
1586e7058c64SPeter Brune    Level: advanced
1587e7058c64SPeter Brune 
1588f1c6b773SPeter Brune .keywords: SNESLineSearch, append, options, prefix, database
1589e7058c64SPeter Brune 
1590e7058c64SPeter Brune .seealso: SNESGetOptionsPrefix()
1591e7058c64SPeter Brune @*/
1592f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[])
1593e7058c64SPeter Brune {
1594e7058c64SPeter Brune   PetscErrorCode ierr;
1595e7058c64SPeter Brune 
1596e7058c64SPeter Brune   PetscFunctionBegin;
1597f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1598e7058c64SPeter Brune   ierr = PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1599e7058c64SPeter Brune   PetscFunctionReturn(0);
1600e7058c64SPeter Brune }
1601e7058c64SPeter Brune 
1602e7058c64SPeter Brune /*@C
1603f1c6b773SPeter Brune    SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all
1604f1c6b773SPeter Brune    SNESLineSearch options in the database.
1605e7058c64SPeter Brune 
1606e7058c64SPeter Brune    Not Collective
1607e7058c64SPeter Brune 
1608e7058c64SPeter Brune    Input Parameter:
1609cd7522eaSPeter Brune .  linesearch - the SNESLineSearch context
1610e7058c64SPeter Brune 
1611e7058c64SPeter Brune    Output Parameter:
1612e7058c64SPeter Brune .  prefix - pointer to the prefix string used
1613e7058c64SPeter Brune 
16148e557f58SPeter Brune    Notes:
16158e557f58SPeter Brune    On the fortran side, the user should pass in a string 'prefix' of
1616e7058c64SPeter Brune    sufficient length to hold the prefix.
1617e7058c64SPeter Brune 
1618e7058c64SPeter Brune    Level: advanced
1619e7058c64SPeter Brune 
1620f1c6b773SPeter Brune .keywords: SNESLineSearch, get, options, prefix, database
1621e7058c64SPeter Brune 
1622e7058c64SPeter Brune .seealso: SNESAppendOptionsPrefix()
1623e7058c64SPeter Brune @*/
1624f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[])
1625e7058c64SPeter Brune {
1626e7058c64SPeter Brune   PetscErrorCode ierr;
1627e7058c64SPeter Brune 
1628e7058c64SPeter Brune   PetscFunctionBegin;
1629f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1630e7058c64SPeter Brune   ierr = PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1631e7058c64SPeter Brune   PetscFunctionReturn(0);
1632e7058c64SPeter Brune }
1633bf7f4e0aSPeter Brune 
16348d359177SBarry Smith /*@C
1635fa0ddf94SBarry Smith    SNESLineSearchSetWorkVecs - Gets work vectors for the line search.
1636f40b411bSPeter Brune 
1637f40b411bSPeter Brune    Input Parameter:
1638f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
1639f40b411bSPeter Brune -  nwork - the number of work vectors
1640f40b411bSPeter Brune 
1641f40b411bSPeter Brune    Level: developer
1642f40b411bSPeter Brune 
1643fa0ddf94SBarry Smith    Developers Note: This is PETSC_EXTERN because it may be used by user written plugin SNES implementations
1644cd7522eaSPeter Brune 
1645f1c6b773SPeter Brune .keywords: SNESLineSearch, work, vector
1646f40b411bSPeter Brune 
1647fa0ddf94SBarry Smith .seealso: SNESSetWorkVecs()
1648f40b411bSPeter Brune @*/
1649fa0ddf94SBarry Smith PetscErrorCode  SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork)
1650bf7f4e0aSPeter Brune {
1651bf7f4e0aSPeter Brune   PetscErrorCode ierr;
1652bf388a1fSBarry Smith 
1653bf7f4e0aSPeter Brune   PetscFunctionBegin;
1654bf7f4e0aSPeter Brune   if (linesearch->vec_sol) {
1655bf7f4e0aSPeter Brune     ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr);
16568d359177SBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
1657bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1658bf7f4e0aSPeter Brune }
1659bf7f4e0aSPeter Brune 
1660f40b411bSPeter Brune /*@
1661422a814eSBarry Smith    SNESLineSearchGetReason - Gets the success/failure status of the last line search application
1662f40b411bSPeter Brune 
1663f40b411bSPeter Brune    Input Parameters:
166478bcb3b5SPeter Brune .  linesearch - linesearch context
1665f40b411bSPeter Brune 
1666f40b411bSPeter Brune    Output Parameters:
1667422a814eSBarry Smith .  result - The success or failure status
1668f40b411bSPeter Brune 
1669cd7522eaSPeter Brune    Notes:
1670c98378a5SBarry Smith    This is typically called after SNESLineSearchApply() in order to determine if the line-search failed
1671cd7522eaSPeter Brune    (and set the SNES convergence accordingly).
1672cd7522eaSPeter Brune 
1673f40b411bSPeter Brune    Level: intermediate
1674f40b411bSPeter Brune 
1675422a814eSBarry Smith .seealso: SNESLineSearchSetReason(), SNESLineSearchReason
1676f40b411bSPeter Brune @*/
1677422a814eSBarry Smith PetscErrorCode  SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result)
1678bf7f4e0aSPeter Brune {
1679bf7f4e0aSPeter Brune   PetscFunctionBegin;
1680f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1681422a814eSBarry Smith   PetscValidPointer(result, 2);
1682422a814eSBarry Smith   *result = linesearch->result;
1683bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1684bf7f4e0aSPeter Brune }
1685bf7f4e0aSPeter Brune 
1686f40b411bSPeter Brune /*@
1687422a814eSBarry Smith    SNESLineSearchSetReason - Sets the success/failure status of the last line search application
1688f40b411bSPeter Brune 
1689f40b411bSPeter Brune    Input Parameters:
169078bcb3b5SPeter Brune +  linesearch - linesearch context
1691422a814eSBarry Smith -  result - The success or failure status
1692f40b411bSPeter Brune 
1693cd7522eaSPeter Brune    Notes:
1694cd7522eaSPeter Brune    This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set
1695cd7522eaSPeter Brune    the success or failure of the line search method.
1696cd7522eaSPeter Brune 
169778bcb3b5SPeter Brune    Level: developer
1698f40b411bSPeter Brune 
1699422a814eSBarry Smith .seealso: SNESLineSearchGetSResult()
1700f40b411bSPeter Brune @*/
1701422a814eSBarry Smith PetscErrorCode  SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result)
17026a388c36SPeter Brune {
17036a388c36SPeter Brune   PetscFunctionBegin;
17045fd66863SKarl Rupp   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1705422a814eSBarry Smith   linesearch->result = result;
17066a388c36SPeter Brune   PetscFunctionReturn(0);
17076a388c36SPeter Brune }
17086a388c36SPeter Brune 
17099bd66eb0SPeter Brune /*@C
1710f1c6b773SPeter Brune    SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.
17119bd66eb0SPeter Brune 
17129bd66eb0SPeter Brune    Input Parameters:
17139bd66eb0SPeter Brune +  snes - nonlinear context obtained from SNESCreate()
17149bd66eb0SPeter Brune .  projectfunc - function for projecting the function to the bounds
17159bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
17169bd66eb0SPeter Brune 
17179bd66eb0SPeter Brune    Logically Collective on SNES
17189bd66eb0SPeter Brune 
17199bd66eb0SPeter Brune    Calling sequence of projectfunc:
17209bd66eb0SPeter Brune .vb
17219bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X)
17229bd66eb0SPeter Brune .ve
17239bd66eb0SPeter Brune 
17249bd66eb0SPeter Brune     Input parameters for projectfunc:
17259bd66eb0SPeter Brune +   snes - nonlinear context
17269bd66eb0SPeter Brune -   X - current solution
17279bd66eb0SPeter Brune 
1728cd7522eaSPeter Brune     Output parameters for projectfunc:
17299bd66eb0SPeter Brune .   X - Projected solution
17309bd66eb0SPeter Brune 
17319bd66eb0SPeter Brune    Calling sequence of normfunc:
17329bd66eb0SPeter Brune .vb
17339bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm)
17349bd66eb0SPeter Brune .ve
17359bd66eb0SPeter Brune 
1736cd7522eaSPeter Brune     Input parameters for normfunc:
17379bd66eb0SPeter Brune +   snes - nonlinear context
17389bd66eb0SPeter Brune .   X - current solution
17399bd66eb0SPeter Brune -   F - current residual
17409bd66eb0SPeter Brune 
1741cd7522eaSPeter Brune     Output parameters for normfunc:
17429bd66eb0SPeter Brune .   fnorm - VI-specific norm of the function
17439bd66eb0SPeter Brune 
1744cd7522eaSPeter Brune     Notes:
1745cd7522eaSPeter Brune     The VI solvers require projection of the solution to the feasible set.  projectfunc should implement this.
1746cd7522eaSPeter Brune 
1747cd7522eaSPeter Brune     The VI solvers require special evaluation of the function norm such that the norm is only calculated
1748cd7522eaSPeter Brune     on the inactive set.  This should be implemented by normfunc.
17499bd66eb0SPeter Brune 
17509bd66eb0SPeter Brune     Level: developer
17519bd66eb0SPeter Brune 
17529bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, set, line search
17539bd66eb0SPeter Brune 
1754f1c6b773SPeter Brune .seealso: SNESLineSearchGetVIFunctions(), SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck()
17559bd66eb0SPeter Brune @*/
1756f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc)
17579bd66eb0SPeter Brune {
17589bd66eb0SPeter Brune   PetscFunctionBegin;
1759f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
17609bd66eb0SPeter Brune   if (projectfunc) linesearch->ops->viproject = projectfunc;
17619bd66eb0SPeter Brune   if (normfunc) linesearch->ops->vinorm = normfunc;
17629bd66eb0SPeter Brune   PetscFunctionReturn(0);
17639bd66eb0SPeter Brune }
17649bd66eb0SPeter Brune 
17659bd66eb0SPeter Brune /*@C
1766f1c6b773SPeter Brune    SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.
17679bd66eb0SPeter Brune 
17689bd66eb0SPeter Brune    Input Parameters:
1769907376e6SBarry Smith .  linesearch - the line search context, obtain with SNESGetLineSearch()
17709bd66eb0SPeter Brune 
17719bd66eb0SPeter Brune    Output Parameters:
17729bd66eb0SPeter Brune +  projectfunc - function for projecting the function to the bounds
17739bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
17749bd66eb0SPeter Brune 
17759bd66eb0SPeter Brune    Logically Collective on SNES
17769bd66eb0SPeter Brune 
17779bd66eb0SPeter Brune     Level: developer
17789bd66eb0SPeter Brune 
17799bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, get, line search
17809bd66eb0SPeter Brune 
1781f1c6b773SPeter Brune .seealso: SNESLineSearchSetVIFunctions(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
17829bd66eb0SPeter Brune @*/
1783f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc)
17849bd66eb0SPeter Brune {
17859bd66eb0SPeter Brune   PetscFunctionBegin;
17869bd66eb0SPeter Brune   if (projectfunc) *projectfunc = linesearch->ops->viproject;
17879bd66eb0SPeter Brune   if (normfunc) *normfunc = linesearch->ops->vinorm;
17889bd66eb0SPeter Brune   PetscFunctionReturn(0);
17899bd66eb0SPeter Brune }
17909bd66eb0SPeter Brune 
1791bf7f4e0aSPeter Brune /*@C
17921c84c290SBarry Smith   SNESLineSearchRegister - See SNESLineSearchRegister()
1793bf7f4e0aSPeter Brune 
1794bf7f4e0aSPeter Brune   Level: advanced
1795bf7f4e0aSPeter Brune @*/
1796bdf89e91SBarry Smith PetscErrorCode  SNESLineSearchRegister(const char sname[],PetscErrorCode (*function)(SNESLineSearch))
1797bf7f4e0aSPeter Brune {
1798bf7f4e0aSPeter Brune   PetscErrorCode ierr;
1799bf7f4e0aSPeter Brune 
1800bf7f4e0aSPeter Brune   PetscFunctionBegin;
1801a240a19fSJed Brown   ierr = PetscFunctionListAdd(&SNESLineSearchList,sname,function);CHKERRQ(ierr);
1802bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1803bf7f4e0aSPeter Brune }
1804