xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision dcf2fd19addf68f961057cf0cf7f6217b6b499d5)
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;
7f1c6b773SPeter Brune PetscLogEvent SNESLineSearch_Apply;
8bf7f4e0aSPeter Brune 
9bf7f4e0aSPeter Brune #undef __FUNCT__
10*dcf2fd19SBarry Smith #define __FUNCT__ "SNESLineSearchMonitorCancel"
11*dcf2fd19SBarry Smith /*@
12*dcf2fd19SBarry Smith    SNESLineSearchMonitorCancel - Clears all the monitor functions for a SNESLineSearch object.
13*dcf2fd19SBarry Smith 
14*dcf2fd19SBarry Smith    Logically Collective on SNESLineSearch
15*dcf2fd19SBarry Smith 
16*dcf2fd19SBarry Smith    Input Parameters:
17*dcf2fd19SBarry Smith .  ls - the SNESLineSearch context
18*dcf2fd19SBarry Smith 
19*dcf2fd19SBarry Smith    Options Database Key:
20*dcf2fd19SBarry Smith .  -snes_linesearch_monitor_cancel - cancels all monitors that have been hardwired
21*dcf2fd19SBarry Smith     into a code by calls to SNESLineSearchMonitorSet(), but does not cancel those
22*dcf2fd19SBarry Smith     set via the options database
23*dcf2fd19SBarry Smith 
24*dcf2fd19SBarry Smith    Notes:
25*dcf2fd19SBarry Smith    There is no way to clear one specific monitor from a SNESLineSearch object.
26*dcf2fd19SBarry Smith 
27*dcf2fd19SBarry Smith    This does not clear the monitor set with SNESLineSearchSetDefaultMonitor() use SNESLineSearchSetDefaultMonitor(ls,NULL) to cancel
28*dcf2fd19SBarry Smith    that one.
29*dcf2fd19SBarry Smith 
30*dcf2fd19SBarry Smith    Level: intermediate
31*dcf2fd19SBarry Smith 
32*dcf2fd19SBarry Smith .keywords: SNESLineSearch, nonlinear, set, monitor
33*dcf2fd19SBarry Smith 
34*dcf2fd19SBarry Smith .seealso: SNESLineSearchMonitorDefault(), SNESLineSearchMonitorSet()
35*dcf2fd19SBarry Smith @*/
36*dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchMonitorCancel(SNESLineSearch ls)
37*dcf2fd19SBarry Smith {
38*dcf2fd19SBarry Smith   PetscErrorCode ierr;
39*dcf2fd19SBarry Smith   PetscInt       i;
40*dcf2fd19SBarry Smith 
41*dcf2fd19SBarry Smith   PetscFunctionBegin;
42*dcf2fd19SBarry Smith   PetscValidHeaderSpecific(ls,SNESLINESEARCH_CLASSID,1);
43*dcf2fd19SBarry Smith   for (i=0; i<ls->numbermonitors; i++) {
44*dcf2fd19SBarry Smith     if (ls->monitordestroy[i]) {
45*dcf2fd19SBarry Smith       ierr = (*ls->monitordestroy[i])(&ls->monitorcontext[i]);CHKERRQ(ierr);
46*dcf2fd19SBarry Smith     }
47*dcf2fd19SBarry Smith   }
48*dcf2fd19SBarry Smith   ls->numbermonitors = 0;
49*dcf2fd19SBarry Smith   PetscFunctionReturn(0);
50*dcf2fd19SBarry Smith }
51*dcf2fd19SBarry Smith 
52*dcf2fd19SBarry Smith #undef __FUNCT__
53*dcf2fd19SBarry Smith #define __FUNCT__ "SNESLineSearchMonitor"
54*dcf2fd19SBarry Smith /*@
55*dcf2fd19SBarry Smith    SNESLineSearchMonitor - runs the user provided monitor routines, if they exist
56*dcf2fd19SBarry Smith 
57*dcf2fd19SBarry Smith    Collective on SNES
58*dcf2fd19SBarry Smith 
59*dcf2fd19SBarry Smith    Input Parameters:
60*dcf2fd19SBarry Smith .  ls - the linesearch object
61*dcf2fd19SBarry Smith 
62*dcf2fd19SBarry Smith    Notes:
63*dcf2fd19SBarry Smith    This routine is called by the SNES implementations.
64*dcf2fd19SBarry Smith    It does not typically need to be called by the user.
65*dcf2fd19SBarry Smith 
66*dcf2fd19SBarry Smith    Level: developer
67*dcf2fd19SBarry Smith 
68*dcf2fd19SBarry Smith .seealso: SNESLineSearchMonitorSet()
69*dcf2fd19SBarry Smith @*/
70*dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchMonitor(SNESLineSearch ls)
71*dcf2fd19SBarry Smith {
72*dcf2fd19SBarry Smith   PetscErrorCode ierr;
73*dcf2fd19SBarry Smith   PetscInt       i,n = ls->numbermonitors;
74*dcf2fd19SBarry Smith 
75*dcf2fd19SBarry Smith   PetscFunctionBegin;
76*dcf2fd19SBarry Smith   for (i=0; i<n; i++) {
77*dcf2fd19SBarry Smith     ierr = (*ls->monitorftns[i])(ls,ls->monitorcontext[i]);CHKERRQ(ierr);
78*dcf2fd19SBarry Smith   }
79*dcf2fd19SBarry Smith   PetscFunctionReturn(0);
80*dcf2fd19SBarry Smith }
81*dcf2fd19SBarry Smith 
82*dcf2fd19SBarry Smith #undef __FUNCT__
83*dcf2fd19SBarry Smith #define __FUNCT__ "SNESLineSearchMonitorSet"
84*dcf2fd19SBarry Smith /*@C
85*dcf2fd19SBarry Smith    SNESLineSearchMonitorSet - Sets an ADDITIONAL function that is to be used at every
86*dcf2fd19SBarry Smith    iteration of the nonlinear solver to display the iteration's
87*dcf2fd19SBarry Smith    progress.
88*dcf2fd19SBarry Smith 
89*dcf2fd19SBarry Smith    Logically Collective on SNESLineSearch
90*dcf2fd19SBarry Smith 
91*dcf2fd19SBarry Smith    Input Parameters:
92*dcf2fd19SBarry Smith +  ls - the SNESLineSearch context
93*dcf2fd19SBarry Smith .  f - the monitor function
94*dcf2fd19SBarry Smith .  mctx - [optional] user-defined context for private data for the
95*dcf2fd19SBarry Smith           monitor routine (use NULL if no context is desired)
96*dcf2fd19SBarry Smith -  monitordestroy - [optional] routine that frees monitor context
97*dcf2fd19SBarry Smith           (may be NULL)
98*dcf2fd19SBarry Smith 
99*dcf2fd19SBarry Smith    Notes:
100*dcf2fd19SBarry Smith    Several different monitoring routines may be set by calling
101*dcf2fd19SBarry Smith    SNESLineSearchMonitorSet() multiple times; all will be called in the
102*dcf2fd19SBarry Smith    order in which they were set.
103*dcf2fd19SBarry Smith 
104*dcf2fd19SBarry Smith    Fortran notes: Only a single monitor function can be set for each SNESLineSearch object
105*dcf2fd19SBarry Smith 
106*dcf2fd19SBarry Smith    Level: intermediate
107*dcf2fd19SBarry Smith 
108*dcf2fd19SBarry Smith .keywords: SNESLineSearch, nonlinear, set, monitor
109*dcf2fd19SBarry Smith 
110*dcf2fd19SBarry Smith .seealso: SNESLineSearchMonitorDefault(), SNESLineSearchMonitorCancel()
111*dcf2fd19SBarry Smith @*/
112*dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchMonitorSet(SNESLineSearch ls,PetscErrorCode (*f)(SNESLineSearch,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
113*dcf2fd19SBarry Smith {
114*dcf2fd19SBarry Smith   PetscFunctionBegin;
115*dcf2fd19SBarry Smith   PetscValidHeaderSpecific(ls,SNESLINESEARCH_CLASSID,1);
116*dcf2fd19SBarry Smith   if (ls->numbermonitors >= MAXSNESLSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
117*dcf2fd19SBarry Smith   ls->monitorftns[ls->numbermonitors]          = f;
118*dcf2fd19SBarry Smith   ls->monitordestroy[ls->numbermonitors]   = monitordestroy;
119*dcf2fd19SBarry Smith   ls->monitorcontext[ls->numbermonitors++] = (void*)mctx;
120*dcf2fd19SBarry Smith   PetscFunctionReturn(0);
121*dcf2fd19SBarry Smith }
122*dcf2fd19SBarry Smith 
123*dcf2fd19SBarry Smith #undef __FUNCT__
124*dcf2fd19SBarry Smith #define __FUNCT__ "SNESLineSearchMonitorSolutionUpdate"
125*dcf2fd19SBarry Smith /*@C
126*dcf2fd19SBarry Smith    SNESLineSearchMonitorSolutionUpdate - Monitors each update a new function value the linesearch tries
127*dcf2fd19SBarry Smith 
128*dcf2fd19SBarry Smith    Collective on SNESLineSearch
129*dcf2fd19SBarry Smith 
130*dcf2fd19SBarry Smith    Input Parameters:
131*dcf2fd19SBarry Smith +  ls - the SNES linesearch object
132*dcf2fd19SBarry Smith -  dummy - the context for the monitor, in this case it is an ASCII PetscViewer
133*dcf2fd19SBarry Smith 
134*dcf2fd19SBarry Smith    Level: intermediate
135*dcf2fd19SBarry Smith 
136*dcf2fd19SBarry Smith .keywords: SNES, nonlinear, default, monitor, norm
137*dcf2fd19SBarry Smith 
138*dcf2fd19SBarry Smith .seealso: SNESMonitorSet(), SNESMonitorSolution()
139*dcf2fd19SBarry Smith @*/
140*dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchMonitorSolutionUpdate(SNESLineSearch ls,void *dummy)
141*dcf2fd19SBarry Smith {
142*dcf2fd19SBarry Smith   PetscErrorCode ierr;
143*dcf2fd19SBarry Smith   PetscViewer    viewer = (PetscViewer) dummy;
144*dcf2fd19SBarry Smith   Vec            Y,W,G;
145*dcf2fd19SBarry Smith 
146*dcf2fd19SBarry Smith   PetscFunctionBegin;
147*dcf2fd19SBarry Smith   ierr = SNESLineSearchGetVecs(ls,NULL,NULL,&Y,&W,&G);CHKERRQ(ierr);
148*dcf2fd19SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"LineSearch attempted update to solution \n");CHKERRQ(ierr);
149*dcf2fd19SBarry Smith   ierr = VecView(Y,viewer);CHKERRQ(ierr);
150*dcf2fd19SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"LineSearch attempted new solution \n");CHKERRQ(ierr);
151*dcf2fd19SBarry Smith   ierr = VecView(W,viewer);CHKERRQ(ierr);
152*dcf2fd19SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"LineSearch attempted updated function value\n");CHKERRQ(ierr);
153*dcf2fd19SBarry Smith   ierr = VecView(G,viewer);CHKERRQ(ierr);
154*dcf2fd19SBarry Smith   PetscFunctionReturn(0);
155*dcf2fd19SBarry Smith }
156*dcf2fd19SBarry Smith 
157*dcf2fd19SBarry Smith #undef __FUNCT__
158f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate"
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 
220bf7f4e0aSPeter Brune #undef __FUNCT__
221f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetUp"
222f40b411bSPeter Brune /*@
22378bcb3b5SPeter Brune    SNESLineSearchSetUp - Prepares the line search for being applied by allocating
22478bcb3b5SPeter Brune    any required vectors.
225f40b411bSPeter Brune 
226cd7522eaSPeter Brune    Collective on SNESLineSearch
227f40b411bSPeter Brune 
228f40b411bSPeter Brune    Input Parameters:
229f40b411bSPeter Brune .  linesearch - The LineSearch instance.
230f40b411bSPeter Brune 
231cd7522eaSPeter Brune    Notes:
232f190f2fcSBarry Smith    For most cases, this needn't be called by users or outside of SNESLineSearchApply().
233cd7522eaSPeter Brune    The only current case where this is called outside of this is for the VI
23478bcb3b5SPeter Brune    solvers, which modify the solution and work vectors before the first call
235cd7522eaSPeter Brune    of SNESLineSearchApply, requiring the SNESLineSearch work vectors to be
236cd7522eaSPeter Brune    allocated upfront.
237cd7522eaSPeter Brune 
23878bcb3b5SPeter Brune    Level: advanced
239f40b411bSPeter Brune 
240f1c6b773SPeter Brune .keywords: SNESLineSearch, SetUp
241f40b411bSPeter Brune 
242f1c6b773SPeter Brune .seealso: SNESLineSearchReset()
243f40b411bSPeter Brune @*/
244f40b411bSPeter Brune 
245bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch)
246bf388a1fSBarry Smith {
247bf7f4e0aSPeter Brune   PetscErrorCode ierr;
248bf388a1fSBarry Smith 
249bf7f4e0aSPeter Brune   PetscFunctionBegin;
250bf7f4e0aSPeter Brune   if (!((PetscObject)linesearch)->type_name) {
2511a4f838cSPeter Brune     ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr);
252bf7f4e0aSPeter Brune   }
253bf7f4e0aSPeter Brune   if (!linesearch->setupcalled) {
2549bd66eb0SPeter Brune     if (!linesearch->vec_sol_new) {
255bf7f4e0aSPeter Brune       ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new);CHKERRQ(ierr);
2569bd66eb0SPeter Brune     }
2579bd66eb0SPeter Brune     if (!linesearch->vec_func_new) {
2589bd66eb0SPeter Brune       ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new);CHKERRQ(ierr);
2599bd66eb0SPeter Brune     }
260bf7f4e0aSPeter Brune     if (linesearch->ops->setup) {
261bf7f4e0aSPeter Brune       ierr = (*linesearch->ops->setup)(linesearch);CHKERRQ(ierr);
262bf7f4e0aSPeter Brune     }
263ed07d7d7SPeter Brune     if (!linesearch->ops->snesfunc) {ierr = SNESLineSearchSetFunction(linesearch,SNESComputeFunction);CHKERRQ(ierr);}
264bf7f4e0aSPeter Brune     linesearch->lambda      = linesearch->damping;
265bf7f4e0aSPeter Brune     linesearch->setupcalled = PETSC_TRUE;
266bf7f4e0aSPeter Brune   }
267bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
268bf7f4e0aSPeter Brune }
269bf7f4e0aSPeter Brune 
270bf7f4e0aSPeter Brune #undef __FUNCT__
271f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchReset"
272f40b411bSPeter Brune 
273f40b411bSPeter Brune /*@
274f190f2fcSBarry Smith    SNESLineSearchReset - Undoes the SNESLineSearchSetUp() and deletes any Vecs or Mats allocated by the line search.
275f40b411bSPeter Brune 
276f1c6b773SPeter Brune    Collective on SNESLineSearch
277f40b411bSPeter Brune 
278f40b411bSPeter Brune    Input Parameters:
279f40b411bSPeter Brune .  linesearch - The LineSearch instance.
280f40b411bSPeter Brune 
281f190f2fcSBarry Smith    Notes: Usually only called by SNESReset()
282f190f2fcSBarry Smith 
283f190f2fcSBarry Smith    Level: developer
284f40b411bSPeter Brune 
285cd7522eaSPeter Brune .keywords: SNESLineSearch, Reset
286f40b411bSPeter Brune 
287f1c6b773SPeter Brune .seealso: SNESLineSearchSetUp()
288f40b411bSPeter Brune @*/
289f40b411bSPeter Brune 
290bf388a1fSBarry Smith PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch)
291bf388a1fSBarry Smith {
292bf7f4e0aSPeter Brune   PetscErrorCode ierr;
293bf388a1fSBarry Smith 
294bf7f4e0aSPeter Brune   PetscFunctionBegin;
295f5af7f23SKarl Rupp   if (linesearch->ops->reset) (*linesearch->ops->reset)(linesearch);
296f5af7f23SKarl Rupp 
297bf7f4e0aSPeter Brune   ierr = VecDestroy(&linesearch->vec_sol_new);CHKERRQ(ierr);
298bf7f4e0aSPeter Brune   ierr = VecDestroy(&linesearch->vec_func_new);CHKERRQ(ierr);
299bf7f4e0aSPeter Brune 
300bf7f4e0aSPeter Brune   ierr = VecDestroyVecs(linesearch->nwork, &linesearch->work);CHKERRQ(ierr);
301f5af7f23SKarl Rupp 
302bf7f4e0aSPeter Brune   linesearch->nwork       = 0;
303bf7f4e0aSPeter Brune   linesearch->setupcalled = PETSC_FALSE;
304bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
305bf7f4e0aSPeter Brune }
306bf7f4e0aSPeter Brune 
307ed07d7d7SPeter Brune #undef __FUNCT__
308ed07d7d7SPeter Brune #define __FUNCT__ "SNESLineSearchSetFunction"
309ed07d7d7SPeter Brune /*@C
310f190f2fcSBarry Smith    SNESLineSearchSetFunction - Sets the function evaluation used by the SNES line search
311ed07d7d7SPeter Brune 
312ed07d7d7SPeter Brune    Input Parameters:
313ed07d7d7SPeter Brune .  linesearch - the SNESLineSearch context
314f190f2fcSBarry Smith +  func       - function evaluation routine
315ed07d7d7SPeter Brune 
316ed07d7d7SPeter Brune    Level: developer
317ed07d7d7SPeter Brune 
318f190f2fcSBarry Smith    Notes: This is used internally by PETSc and not called by users
319f190f2fcSBarry Smith 
320ed07d7d7SPeter Brune .keywords: get, linesearch, pre-check
321ed07d7d7SPeter Brune 
322f190f2fcSBarry Smith .seealso: SNESSetFunction()
323ed07d7d7SPeter Brune @*/
324ed07d7d7SPeter Brune PetscErrorCode  SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES,Vec,Vec))
325ed07d7d7SPeter Brune {
326ed07d7d7SPeter Brune   PetscFunctionBegin;
327ed07d7d7SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
328ed07d7d7SPeter Brune   linesearch->ops->snesfunc = func;
329ed07d7d7SPeter Brune   PetscFunctionReturn(0);
330ed07d7d7SPeter Brune }
331ed07d7d7SPeter Brune 
332ed07d7d7SPeter Brune 
3336b2b7091SBarry Smith /*MC
334f190f2fcSBarry Smith     SNESLineSearchPreCheckFunction - form of function passed to check the search direction before line search is called
3356b2b7091SBarry Smith 
3366b2b7091SBarry Smith      Synopsis:
337aaa7dc30SBarry Smith      #include <petscsnes.h>
3386b2b7091SBarry Smith      SNESLineSearchPreCheckFunction(SNESLineSearch snes,Vec x,Vec y, PetscBool *changed);
3396b2b7091SBarry Smith 
3406b2b7091SBarry Smith        Input Parameters:
3416b2b7091SBarry Smith +      x - solution vector
3426b2b7091SBarry Smith .      y - search direction vector
3436b2b7091SBarry Smith -      changed - flag to indicate the precheck changed x or y.
3446b2b7091SBarry Smith 
345f190f2fcSBarry Smith      Note: This is NOTE a PETSc function, rather it documents the calling sequence of functions passed to SNESLineSearchSetPreCheck()
346f190f2fcSBarry Smith            and SNESLineSearchGetPreCheck()
347f190f2fcSBarry Smith 
348878cb397SSatish Balay    Level: advanced
349878cb397SSatish Balay 
350f190f2fcSBarry Smith .seealso:   SNESLineSearchSetPreCheck(), SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck()
3516b2b7091SBarry Smith M*/
3526b2b7091SBarry Smith 
35386d74e61SPeter Brune #undef __FUNCT__
354f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetPreCheck"
35586d74e61SPeter Brune /*@C
356f190f2fcSBarry Smith    SNESLineSearchSetPreCheck - Sets a user function that is called after the initial search direction has been computed but
357df3898eeSBarry Smith          before the line search routine has been applied. Allows the user to adjust the result of (usually a linear solve) that
358f190f2fcSBarry Smith          determined the search direction.
35986d74e61SPeter Brune 
360f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
36186d74e61SPeter Brune 
36286d74e61SPeter Brune    Input Parameters:
363f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
364f190f2fcSBarry Smith .  func - [optional] function evaluation routine, see SNESLineSearchPreCheckFunction for the calling sequence
365f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
36686d74e61SPeter Brune 
36786d74e61SPeter Brune    Level: intermediate
36886d74e61SPeter Brune 
36986d74e61SPeter Brune .keywords: set, linesearch, pre-check
37086d74e61SPeter Brune 
371f190f2fcSBarry Smith .seealso: SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
37286d74e61SPeter Brune @*/
373f190f2fcSBarry Smith PetscErrorCode  SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void *ctx)
37486d74e61SPeter Brune {
3759bd66eb0SPeter Brune   PetscFunctionBegin;
376f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
377f190f2fcSBarry Smith   if (func) linesearch->ops->precheck = func;
37886d74e61SPeter Brune   if (ctx) linesearch->precheckctx = ctx;
37986d74e61SPeter Brune   PetscFunctionReturn(0);
38086d74e61SPeter Brune }
38186d74e61SPeter Brune 
38286d74e61SPeter Brune #undef __FUNCT__
383f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetPreCheck"
38486d74e61SPeter Brune /*@C
38553deb899SBarry Smith    SNESLineSearchGetPreCheck - Gets the pre-check function for the line search routine.
38686d74e61SPeter Brune 
38786d74e61SPeter Brune    Input Parameters:
388f1c6b773SPeter Brune .  linesearch - the SNESLineSearch context
38986d74e61SPeter Brune 
39086d74e61SPeter Brune    Output Parameters:
391f190f2fcSBarry Smith +  func       - [optional] function evaluation routine, see SNESLineSearchPreCheckFunction for calling sequence
392f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
39386d74e61SPeter Brune 
39486d74e61SPeter Brune    Level: intermediate
39586d74e61SPeter Brune 
39686d74e61SPeter Brune .keywords: get, linesearch, pre-check
39786d74e61SPeter Brune 
398f1c6b773SPeter Brune .seealso: SNESLineSearchGetPostCheck(), SNESLineSearchSetPreCheck()
39986d74e61SPeter Brune @*/
400f190f2fcSBarry Smith PetscErrorCode  SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void **ctx)
40186d74e61SPeter Brune {
4029bd66eb0SPeter Brune   PetscFunctionBegin;
403f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
404f190f2fcSBarry Smith   if (func) *func = linesearch->ops->precheck;
40586d74e61SPeter Brune   if (ctx) *ctx = linesearch->precheckctx;
40686d74e61SPeter Brune   PetscFunctionReturn(0);
40786d74e61SPeter Brune }
40886d74e61SPeter Brune 
4096b2b7091SBarry Smith /*MC
410f190f2fcSBarry Smith     SNESLineSearchPostCheckFunction - form of function that is called after line search is complete
4116b2b7091SBarry Smith 
4126b2b7091SBarry Smith      Synopsis:
413aaa7dc30SBarry Smith      #include <petscsnes.h>
4146b2b7091SBarry Smith      SNESLineSearchPostheckFunction(SNESLineSearch linesearch,Vec x,Vec y,  Vec w, *changed_y, PetscBool *changed_w);
4156b2b7091SBarry Smith 
4166b2b7091SBarry Smith      Input Parameters:
4176b2b7091SBarry Smith +      x - old solution vector
4186b2b7091SBarry Smith .      y - search direction vector
4196b2b7091SBarry Smith .      w - new solution vector
4206b2b7091SBarry Smith .      changed_y - indicates that the line search changed y
4216b2b7091SBarry Smith -      changed_w - indicates that the line search changed w
4226b2b7091SBarry Smith 
423f190f2fcSBarry Smith      Note: This is NOTE a PETSc function, rather it documents the calling sequence of functions passed to SNESLineSearchSetPostCheck()
424f190f2fcSBarry Smith            and SNESLineSearchGetPostCheck()
425f190f2fcSBarry Smith 
426878cb397SSatish Balay    Level: advanced
4276b2b7091SBarry Smith 
428f190f2fcSBarry Smith .seealso:   SNESLineSearchSetPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPreCheck(), SNESLineSearchGetPostCheck()
4296b2b7091SBarry Smith M*/
4306b2b7091SBarry Smith 
43186d74e61SPeter Brune #undef __FUNCT__
432f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetPostCheck"
43386d74e61SPeter Brune /*@C
434f190f2fcSBarry Smith    SNESLineSearchSetPostCheck - Sets a user function that is called after the line search has been applied to determine the step
435f190f2fcSBarry Smith        direction and length. Allows the user a chance to change or override the decision of the line search routine
43686d74e61SPeter Brune 
437f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
43886d74e61SPeter Brune 
43986d74e61SPeter Brune    Input Parameters:
440f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
441f190f2fcSBarry Smith .  func - [optional] function evaluation routine, see SNESLineSearchPostCheckFunction for the calling sequence
442f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
44386d74e61SPeter Brune 
44486d74e61SPeter Brune    Level: intermediate
44586d74e61SPeter Brune 
44686d74e61SPeter Brune .keywords: set, linesearch, post-check
44786d74e61SPeter Brune 
448f1c6b773SPeter Brune .seealso: SNESLineSearchSetPreCheck()
44986d74e61SPeter Brune @*/
450f190f2fcSBarry Smith PetscErrorCode  SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void *ctx)
45186d74e61SPeter Brune {
45286d74e61SPeter Brune   PetscFunctionBegin;
453f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
454f190f2fcSBarry Smith   if (func) linesearch->ops->postcheck = func;
45586d74e61SPeter Brune   if (ctx) linesearch->postcheckctx = ctx;
45686d74e61SPeter Brune   PetscFunctionReturn(0);
45786d74e61SPeter Brune }
45886d74e61SPeter Brune 
45986d74e61SPeter Brune #undef __FUNCT__
460f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetPostCheck"
46186d74e61SPeter Brune /*@C
462f1c6b773SPeter Brune    SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine.
46386d74e61SPeter Brune 
46486d74e61SPeter Brune    Input Parameters:
465f1c6b773SPeter Brune .  linesearch - the SNESLineSearch context
46686d74e61SPeter Brune 
46786d74e61SPeter Brune    Output Parameters:
468f190f2fcSBarry Smith +  func - [optional] function evaluation routine, see for the calling sequence SNESLineSearchPostCheckFunction
469f190f2fcSBarry Smith -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
47086d74e61SPeter Brune 
47186d74e61SPeter Brune    Level: intermediate
47286d74e61SPeter Brune 
47386d74e61SPeter Brune .keywords: get, linesearch, post-check
47486d74e61SPeter Brune 
475f1c6b773SPeter Brune .seealso: SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck()
47686d74e61SPeter Brune @*/
477f190f2fcSBarry Smith PetscErrorCode  SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void **ctx)
47886d74e61SPeter Brune {
4799bd66eb0SPeter Brune   PetscFunctionBegin;
480f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
481f190f2fcSBarry Smith   if (func) *func = linesearch->ops->postcheck;
48286d74e61SPeter Brune   if (ctx) *ctx = linesearch->postcheckctx;
48386d74e61SPeter Brune   PetscFunctionReturn(0);
48486d74e61SPeter Brune }
48586d74e61SPeter Brune 
486bf7f4e0aSPeter Brune #undef __FUNCT__
487f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPreCheck"
488f40b411bSPeter Brune /*@
489f1c6b773SPeter Brune    SNESLineSearchPreCheck - Prepares the line search for being applied.
490f40b411bSPeter Brune 
491cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
492f40b411bSPeter Brune 
493f40b411bSPeter Brune    Input Parameters:
4947b1df9c1SPeter Brune +  linesearch - The linesearch instance.
4957b1df9c1SPeter Brune .  X - The current solution
4967b1df9c1SPeter Brune -  Y - The step direction
497f40b411bSPeter Brune 
498f40b411bSPeter Brune    Output Parameters:
4998e557f58SPeter Brune .  changed - Indicator that the precheck routine has changed anything
500f40b411bSPeter Brune 
501f190f2fcSBarry Smith    Level: developer
502f40b411bSPeter Brune 
503f1c6b773SPeter Brune .keywords: SNESLineSearch, Create
504f40b411bSPeter Brune 
505f1c6b773SPeter Brune .seealso: SNESLineSearchPostCheck()
506f40b411bSPeter Brune @*/
5077b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed)
508bf7f4e0aSPeter Brune {
509bf7f4e0aSPeter Brune   PetscErrorCode ierr;
5105fd66863SKarl Rupp 
511bf7f4e0aSPeter Brune   PetscFunctionBegin;
512bf7f4e0aSPeter Brune   *changed = PETSC_FALSE;
5136b2b7091SBarry Smith   if (linesearch->ops->precheck) {
5146b2b7091SBarry Smith     ierr = (*linesearch->ops->precheck)(linesearch, X, Y, changed, linesearch->precheckctx);CHKERRQ(ierr);
51538bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed,4);
516bf7f4e0aSPeter Brune   }
517bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
518bf7f4e0aSPeter Brune }
519bf7f4e0aSPeter Brune 
520bf7f4e0aSPeter Brune #undef __FUNCT__
521f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPostCheck"
522f40b411bSPeter Brune /*@
523f1c6b773SPeter Brune    SNESLineSearchPostCheck - Prepares the line search for being applied.
524f40b411bSPeter Brune 
525cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
526f40b411bSPeter Brune 
527f40b411bSPeter Brune    Input Parameters:
5287b1df9c1SPeter Brune +  linesearch - The linesearch context
5297b1df9c1SPeter Brune .  X - The last solution
5307b1df9c1SPeter Brune .  Y - The step direction
5317b1df9c1SPeter Brune -  W - The updated solution, W = X + lambda*Y for some lambda
532f40b411bSPeter Brune 
533f40b411bSPeter Brune    Output Parameters:
53478bcb3b5SPeter Brune +  changed_Y - Indicator if the direction Y has been changed.
53578bcb3b5SPeter Brune -  changed_W - Indicator if the new candidate solution W has been changed.
536f40b411bSPeter Brune 
537f190f2fcSBarry Smith    Level: developer
538f40b411bSPeter Brune 
539f1c6b773SPeter Brune .keywords: SNESLineSearch, Create
540f40b411bSPeter Brune 
541f1c6b773SPeter Brune .seealso: SNESLineSearchPreCheck()
542f40b411bSPeter Brune @*/
5437b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch,Vec X,Vec Y,Vec W,PetscBool *changed_Y,PetscBool *changed_W)
544bf7f4e0aSPeter Brune {
545bf7f4e0aSPeter Brune   PetscErrorCode ierr;
546bf388a1fSBarry Smith 
547bf7f4e0aSPeter Brune   PetscFunctionBegin;
548bf7f4e0aSPeter Brune   *changed_Y = PETSC_FALSE;
549bf7f4e0aSPeter Brune   *changed_W = PETSC_FALSE;
5506b2b7091SBarry Smith   if (linesearch->ops->postcheck) {
5516b2b7091SBarry Smith     ierr = (*linesearch->ops->postcheck)(linesearch,X,Y,W,changed_Y,changed_W,linesearch->postcheckctx);CHKERRQ(ierr);
55238bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed_Y,5);
55338bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed_W,6);
55486d74e61SPeter Brune   }
55586d74e61SPeter Brune   PetscFunctionReturn(0);
55686d74e61SPeter Brune }
55786d74e61SPeter Brune 
55886d74e61SPeter Brune #undef __FUNCT__
559f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPreCheckPicard"
56086d74e61SPeter Brune /*@C
56186d74e61SPeter Brune    SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration
56286d74e61SPeter Brune 
563cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
56486d74e61SPeter Brune 
56586d74e61SPeter Brune    Input Arguments:
56686d74e61SPeter Brune +  linesearch - linesearch context
56786d74e61SPeter Brune .  X - base state for this step
56886d74e61SPeter Brune .  Y - initial correction
569907376e6SBarry Smith -  ctx - context for this function
57086d74e61SPeter Brune 
57186d74e61SPeter Brune    Output Arguments:
57286d74e61SPeter Brune +  Y - correction, possibly modified
57386d74e61SPeter Brune -  changed - flag indicating that Y was modified
57486d74e61SPeter Brune 
57586d74e61SPeter Brune    Options Database Key:
576cd7522eaSPeter Brune +  -snes_linesearch_precheck_picard - activate this routine
577cd7522eaSPeter Brune -  -snes_linesearch_precheck_picard_angle - angle
57886d74e61SPeter Brune 
57986d74e61SPeter Brune    Level: advanced
58086d74e61SPeter Brune 
58186d74e61SPeter Brune    Notes:
58286d74e61SPeter Brune    This function should be passed to SNESLineSearchSetPreCheck()
58386d74e61SPeter Brune 
58486d74e61SPeter Brune    The justification for this method involves the linear convergence of a Picard iteration
58586d74e61SPeter Brune    so the Picard linearization should be provided in place of the "Jacobian". This correction
58686d74e61SPeter Brune    is generally not useful when using a Newton linearization.
58786d74e61SPeter Brune 
58886d74e61SPeter Brune    Reference:
58986d74e61SPeter Brune    Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology.
59086d74e61SPeter Brune 
59186d74e61SPeter Brune .seealso: SNESLineSearchSetPreCheck()
59286d74e61SPeter Brune @*/
593f1c6b773SPeter Brune PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed,void *ctx)
59486d74e61SPeter Brune {
59586d74e61SPeter Brune   PetscErrorCode ierr;
59686d74e61SPeter Brune   PetscReal      angle = *(PetscReal*)linesearch->precheckctx;
59786d74e61SPeter Brune   Vec            Ylast;
59886d74e61SPeter Brune   PetscScalar    dot;
59986d74e61SPeter Brune   PetscInt       iter;
60086d74e61SPeter Brune   PetscReal      ynorm,ylastnorm,theta,angle_radians;
60186d74e61SPeter Brune   SNES           snes;
60286d74e61SPeter Brune 
60386d74e61SPeter Brune   PetscFunctionBegin;
604f1c6b773SPeter Brune   ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
60586d74e61SPeter Brune   ierr = PetscObjectQuery((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject*)&Ylast);CHKERRQ(ierr);
60686d74e61SPeter Brune   if (!Ylast) {
60786d74e61SPeter Brune     ierr = VecDuplicate(Y,&Ylast);CHKERRQ(ierr);
60886d74e61SPeter Brune     ierr = PetscObjectCompose((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject)Ylast);CHKERRQ(ierr);
60986d74e61SPeter Brune     ierr = PetscObjectDereference((PetscObject)Ylast);CHKERRQ(ierr);
61086d74e61SPeter Brune   }
61186d74e61SPeter Brune   ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr);
61286d74e61SPeter Brune   if (iter < 2) {
61386d74e61SPeter Brune     ierr     = VecCopy(Y,Ylast);CHKERRQ(ierr);
61486d74e61SPeter Brune     *changed = PETSC_FALSE;
61586d74e61SPeter Brune     PetscFunctionReturn(0);
61686d74e61SPeter Brune   }
61786d74e61SPeter Brune 
61886d74e61SPeter Brune   ierr = VecDot(Y,Ylast,&dot);CHKERRQ(ierr);
61986d74e61SPeter Brune   ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr);
62086d74e61SPeter Brune   ierr = VecNorm(Ylast,NORM_2,&ylastnorm);CHKERRQ(ierr);
62186d74e61SPeter Brune   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
622255453a1SBarry Smith   theta         = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0));
62386d74e61SPeter Brune   angle_radians = angle * PETSC_PI / 180.;
62486d74e61SPeter Brune   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
62586d74e61SPeter Brune     /* Modify the step Y */
62686d74e61SPeter Brune     PetscReal alpha,ydiffnorm;
62786d74e61SPeter Brune     ierr  = VecAXPY(Ylast,-1.0,Y);CHKERRQ(ierr);
62886d74e61SPeter Brune     ierr  = VecNorm(Ylast,NORM_2,&ydiffnorm);CHKERRQ(ierr);
62986d74e61SPeter Brune     alpha = ylastnorm / ydiffnorm;
63086d74e61SPeter Brune     ierr  = VecCopy(Y,Ylast);CHKERRQ(ierr);
63186d74e61SPeter Brune     ierr  = VecScale(Y,alpha);CHKERRQ(ierr);
632c69d1a72SBarry 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);
63386d74e61SPeter Brune   } else {
634c69d1a72SBarry 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);
63586d74e61SPeter Brune     ierr     = VecCopy(Y,Ylast);CHKERRQ(ierr);
63686d74e61SPeter Brune     *changed = PETSC_FALSE;
637bf7f4e0aSPeter Brune   }
638bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
639bf7f4e0aSPeter Brune }
640bf7f4e0aSPeter Brune 
641bf7f4e0aSPeter Brune #undef __FUNCT__
642f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchApply"
643f40b411bSPeter Brune /*@
644cd7522eaSPeter Brune    SNESLineSearchApply - Computes the line-search update.
645f40b411bSPeter Brune 
646f1c6b773SPeter Brune    Collective on SNESLineSearch
647f40b411bSPeter Brune 
648f40b411bSPeter Brune    Input Parameters:
6498e557f58SPeter Brune +  linesearch - The linesearch context
6508e557f58SPeter Brune .  X - The current solution
6518e557f58SPeter Brune .  F - The current function
6528e557f58SPeter Brune .  fnorm - The current norm
6538e557f58SPeter Brune -  Y - The search direction
654f40b411bSPeter Brune 
655f40b411bSPeter Brune    Output Parameters:
6568e557f58SPeter Brune +  X - The new solution
6578e557f58SPeter Brune .  F - The new function
6588e557f58SPeter Brune -  fnorm - The new function norm
659f40b411bSPeter Brune 
660cd7522eaSPeter Brune    Options Database Keys:
661d4c6564cSPatrick Farrell + -snes_linesearch_type - basic, bt, l2, cp, nleqerr, shell
662*dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches
6633c7d6663SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter
664cd7522eaSPeter Brune . -snes_linesearch_norms   - Turn on/off the linesearch norms
6653c7d6663SPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess
6663c7d6663SPeter Brune - -snes_linesearch_max_it - The number of iterations for iterative line searches
667cd7522eaSPeter Brune 
668cd7522eaSPeter Brune    Notes:
669cd7522eaSPeter Brune    This is typically called from within a SNESSolve() implementation in order to
670cd7522eaSPeter Brune    help with convergence of the nonlinear method.  Various SNES types use line searches
671cd7522eaSPeter Brune    in different ways, but the overarching theme is that a line search is used to determine
672cd7522eaSPeter Brune    an optimal damping parameter of a step at each iteration of the method.  Each
673cd7522eaSPeter Brune    application of the line search may invoke SNESComputeFunction several times, and
674cd7522eaSPeter Brune    therefore may be fairly expensive.
675cd7522eaSPeter Brune 
676f40b411bSPeter Brune    Level: Intermediate
677f40b411bSPeter Brune 
678f1c6b773SPeter Brune .keywords: SNESLineSearch, Create
679f40b411bSPeter Brune 
680cd7522eaSPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchPreCheck(), SNESLineSearchPostCheck(), SNESSolve(), SNESComputeFunction()
681f40b411bSPeter Brune @*/
682bf388a1fSBarry Smith PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y)
683bf388a1fSBarry Smith {
684bf7f4e0aSPeter Brune   PetscErrorCode ierr;
685bf7f4e0aSPeter Brune 
686bf388a1fSBarry Smith   PetscFunctionBegin;
687f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
688bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
689bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
690bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
691bf7f4e0aSPeter Brune 
692422a814eSBarry Smith   linesearch->result = SNES_LINESEARCH_SUCCEEDED;
693bf7f4e0aSPeter Brune 
694bf7f4e0aSPeter Brune   linesearch->vec_sol    = X;
695bf7f4e0aSPeter Brune   linesearch->vec_update = Y;
696bf7f4e0aSPeter Brune   linesearch->vec_func   = F;
697bf7f4e0aSPeter Brune 
698f1c6b773SPeter Brune   ierr = SNESLineSearchSetUp(linesearch);CHKERRQ(ierr);
699bf7f4e0aSPeter Brune 
700f5af7f23SKarl Rupp   if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */
701bf7f4e0aSPeter Brune 
702f5af7f23SKarl Rupp   if (fnorm) linesearch->fnorm = *fnorm;
703f5af7f23SKarl Rupp   else {
704bf7f4e0aSPeter Brune     ierr = VecNorm(F, NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
705bf7f4e0aSPeter Brune   }
706bf7f4e0aSPeter Brune 
707f1c6b773SPeter Brune   ierr = PetscLogEventBegin(SNESLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
708bf7f4e0aSPeter Brune 
709bf7f4e0aSPeter Brune   ierr = (*linesearch->ops->apply)(linesearch);CHKERRQ(ierr);
710bf7f4e0aSPeter Brune 
711f1c6b773SPeter Brune   ierr = PetscLogEventEnd(SNESLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
712bf7f4e0aSPeter Brune 
713f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
714bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
715bf7f4e0aSPeter Brune }
716bf7f4e0aSPeter Brune 
717bf7f4e0aSPeter Brune #undef __FUNCT__
718f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchDestroy"
719f40b411bSPeter Brune /*@
720f1c6b773SPeter Brune    SNESLineSearchDestroy - Destroys the line search instance.
721f40b411bSPeter Brune 
722f1c6b773SPeter Brune    Collective on SNESLineSearch
723f40b411bSPeter Brune 
724f40b411bSPeter Brune    Input Parameters:
7258e557f58SPeter Brune .  linesearch - The linesearch context
726f40b411bSPeter Brune 
727f40b411bSPeter Brune    Level: Intermediate
728f40b411bSPeter Brune 
72978bcb3b5SPeter Brune .keywords: SNESLineSearch, Destroy
730f40b411bSPeter Brune 
731cd7522eaSPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchReset(), SNESDestroy()
732f40b411bSPeter Brune @*/
733bf388a1fSBarry Smith PetscErrorCode SNESLineSearchDestroy(SNESLineSearch * linesearch)
734bf388a1fSBarry Smith {
735bf7f4e0aSPeter Brune   PetscErrorCode ierr;
736bf388a1fSBarry Smith 
737bf7f4e0aSPeter Brune   PetscFunctionBegin;
738bf7f4e0aSPeter Brune   if (!*linesearch) PetscFunctionReturn(0);
739f1c6b773SPeter Brune   PetscValidHeaderSpecific((*linesearch),SNESLINESEARCH_CLASSID,1);
740bf7f4e0aSPeter Brune   if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = 0; PetscFunctionReturn(0);}
741e04113cfSBarry Smith   ierr = PetscObjectSAWsViewOff((PetscObject)*linesearch);CHKERRQ(ierr);
74222d28d08SBarry Smith   ierr = SNESLineSearchReset(*linesearch);CHKERRQ(ierr);
743f5af7f23SKarl Rupp   if ((*linesearch)->ops->destroy) (*linesearch)->ops->destroy(*linesearch);
744bf7f4e0aSPeter Brune   ierr = PetscViewerDestroy(&(*linesearch)->monitor);CHKERRQ(ierr);
745*dcf2fd19SBarry Smith   ierr = SNESLineSearchMonitorCancel((*linesearch));CHKERRQ(ierr);
746e7058c64SPeter Brune   ierr = PetscHeaderDestroy(linesearch);CHKERRQ(ierr);
747bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
748bf7f4e0aSPeter Brune }
749bf7f4e0aSPeter Brune 
750bf7f4e0aSPeter Brune #undef __FUNCT__
751*dcf2fd19SBarry Smith #define __FUNCT__ "SNESLineSearchSetDefaultMonitor"
752f40b411bSPeter Brune /*@
753*dcf2fd19SBarry Smith    SNESLineSearchSetDefaultMonitor - Turns on/off printing useful information and debugging output about the line search.
754bf7f4e0aSPeter Brune 
755bf7f4e0aSPeter Brune    Input Parameters:
756*dcf2fd19SBarry Smith +  linesearch - the linesearch object
757*dcf2fd19SBarry Smith -  viewer - an ASCII PetscViewer or NULL to turn off monitor
758bf7f4e0aSPeter Brune 
759*dcf2fd19SBarry Smith    Logically Collective on SNESLineSearch
760bf7f4e0aSPeter Brune 
761bf7f4e0aSPeter Brune    Options Database:
762*dcf2fd19SBarry Smith .   -snes_linesearch_monitor [:filename] - enables the monitor
763bf7f4e0aSPeter Brune 
764bf7f4e0aSPeter Brune    Level: intermediate
765bf7f4e0aSPeter Brune 
766*dcf2fd19SBarry Smith    Developer Note: This monitor is implemented differently than the other SNESLineSearchMonitors that are set with
767*dcf2fd19SBarry Smith      SNESLineSearchMonitorSet() since it is called in many locations of the line search routines.
768*dcf2fd19SBarry Smith 
769*dcf2fd19SBarry Smith .seealso: SNESLineSearchGetDefaultMonitor(), PetscViewer
770bf7f4e0aSPeter Brune @*/
771*dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchSetDefaultMonitor(SNESLineSearch linesearch, PetscViewer viewer)
772bf7f4e0aSPeter Brune {
773bf7f4e0aSPeter Brune   PetscErrorCode ierr;
774bf388a1fSBarry Smith 
775bf7f4e0aSPeter Brune   PetscFunctionBegin;
776*dcf2fd19SBarry Smith   if (viewer) {ierr = PetscObjectReference((PetscObject)viewer);CHKERRQ(ierr);}
777bf7f4e0aSPeter Brune   ierr = PetscViewerDestroy(&linesearch->monitor);CHKERRQ(ierr);
778*dcf2fd19SBarry Smith   linesearch->monitor = viewer;
779bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
780bf7f4e0aSPeter Brune }
781bf7f4e0aSPeter Brune 
782bf7f4e0aSPeter Brune #undef __FUNCT__
783*dcf2fd19SBarry Smith #define __FUNCT__ "SNESLineSearchGetDefaultMonitor"
784f40b411bSPeter Brune /*@
785*dcf2fd19SBarry Smith    SNESLineSearchGetDefaultMonitor - Gets the PetscViewer instance for the line search monitor.
7866a388c36SPeter Brune 
787f190f2fcSBarry Smith    Input Parameter:
7888e557f58SPeter Brune .  linesearch - linesearch context
789f40b411bSPeter Brune 
790f190f2fcSBarry Smith    Output Parameter:
7918e557f58SPeter Brune .  monitor - monitor context
792f40b411bSPeter Brune 
793f40b411bSPeter Brune    Logically Collective on SNES
794f40b411bSPeter Brune 
795f40b411bSPeter Brune    Options Database Keys:
7968e557f58SPeter Brune .   -snes_linesearch_monitor - enables the monitor
797f40b411bSPeter Brune 
798f40b411bSPeter Brune    Level: intermediate
799f40b411bSPeter Brune 
800*dcf2fd19SBarry Smith .seealso: SNESLineSearchSetDefaultMonitor(), PetscViewer
801f40b411bSPeter Brune @*/
802*dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchGetDefaultMonitor(SNESLineSearch linesearch, PetscViewer *monitor)
8036a388c36SPeter Brune {
8046a388c36SPeter Brune   PetscFunctionBegin;
805f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
8066a388c36SPeter Brune   if (monitor) {
8076a388c36SPeter Brune     PetscValidPointer(monitor, 2);
8086a388c36SPeter Brune     *monitor = linesearch->monitor;
8096a388c36SPeter Brune   }
8106a388c36SPeter Brune   PetscFunctionReturn(0);
8116a388c36SPeter Brune }
8126a388c36SPeter Brune 
8136a388c36SPeter Brune #undef __FUNCT__
814*dcf2fd19SBarry Smith #define __FUNCT__ "SNESLineSearchMonitorSetFromOptions"
815*dcf2fd19SBarry Smith /*@C
816*dcf2fd19SBarry Smith    SNESLineSearchMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
817*dcf2fd19SBarry Smith 
818*dcf2fd19SBarry Smith    Collective on SNESLineSearch
819*dcf2fd19SBarry Smith 
820*dcf2fd19SBarry Smith    Input Parameters:
821*dcf2fd19SBarry Smith +  ls - LineSearch object you wish to monitor
822*dcf2fd19SBarry Smith .  name - the monitor type one is seeking
823*dcf2fd19SBarry Smith .  help - message indicating what monitoring is done
824*dcf2fd19SBarry Smith .  manual - manual page for the monitor
825*dcf2fd19SBarry Smith .  monitor - the monitor function
826*dcf2fd19SBarry 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
827*dcf2fd19SBarry Smith 
828*dcf2fd19SBarry Smith    Level: developer
829*dcf2fd19SBarry Smith 
830*dcf2fd19SBarry Smith .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
831*dcf2fd19SBarry Smith           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
832*dcf2fd19SBarry Smith           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
833*dcf2fd19SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
834*dcf2fd19SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
835*dcf2fd19SBarry Smith           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
836*dcf2fd19SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
837*dcf2fd19SBarry Smith @*/
838*dcf2fd19SBarry Smith PetscErrorCode  SNESLineSearchMonitorSetFromOptions(SNESLineSearch ls,const char name[],const char help[], const char manual[],PetscErrorCode (*monitor)(SNESLineSearch,void*),PetscErrorCode (*monitorsetup)(SNESLineSearch,PetscViewer))
839*dcf2fd19SBarry Smith {
840*dcf2fd19SBarry Smith   PetscErrorCode    ierr;
841*dcf2fd19SBarry Smith   PetscViewer       viewer;
842*dcf2fd19SBarry Smith   PetscViewerFormat format;
843*dcf2fd19SBarry Smith   PetscBool         flg;
844*dcf2fd19SBarry Smith 
845*dcf2fd19SBarry Smith   PetscFunctionBegin;
846*dcf2fd19SBarry Smith   ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)ls),((PetscObject)ls)->prefix,name,&viewer,&format,&flg);CHKERRQ(ierr);
847*dcf2fd19SBarry Smith   if (flg) {
848*dcf2fd19SBarry Smith     ierr = PetscViewerSetFormat(viewer,format);CHKERRQ(ierr);
849*dcf2fd19SBarry Smith     if (monitorsetup) {
850*dcf2fd19SBarry Smith       ierr = (*monitorsetup)(ls,viewer);CHKERRQ(ierr);
851*dcf2fd19SBarry Smith     }
852*dcf2fd19SBarry Smith     ierr = SNESLineSearchMonitorSet(ls,monitor,viewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
853*dcf2fd19SBarry Smith   }
854*dcf2fd19SBarry Smith   PetscFunctionReturn(0);
855*dcf2fd19SBarry Smith }
856*dcf2fd19SBarry Smith 
857*dcf2fd19SBarry Smith #undef __FUNCT__
858f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetFromOptions"
859f40b411bSPeter Brune /*@
860f1c6b773SPeter Brune    SNESLineSearchSetFromOptions - Sets options for the line search
861f40b411bSPeter Brune 
862f40b411bSPeter Brune    Input Parameters:
8638e557f58SPeter Brune .  linesearch - linesearch context
864f40b411bSPeter Brune 
865f40b411bSPeter Brune    Options Database Keys:
866d4c6564cSPatrick Farrell + -snes_linesearch_type <type> - basic, bt, l2, cp, nleqerr, shell
8673c7d6663SPeter Brune . -snes_linesearch_order <order> - 1, 2, 3.  Most types only support certain orders (bt supports 2 or 3)
8685a9b6599SPeter Brune . -snes_linesearch_norms   - Turn on/off the linesearch norms for the basic linesearch type
86971eef1aeSPeter Brune . -snes_linesearch_minlambda - The minimum step length
8701a9b3a06SPeter Brune . -snes_linesearch_maxstep - The maximum step size
8711a9b3a06SPeter Brune . -snes_linesearch_rtol - Relative tolerance for iterative line searches
8721a9b3a06SPeter Brune . -snes_linesearch_atol - Absolute tolerance for iterative line searches
8731a9b3a06SPeter Brune . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches
8741a9b3a06SPeter Brune . -snes_linesearch_max_it - The number of iterations for iterative line searches
875*dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches
876*dcf2fd19SBarry Smith . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine
8778e557f58SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter
878cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess.
8791a9b3a06SPeter Brune . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method
8801a9b3a06SPeter Brune - -snes_linesearch_precheck_picard_angle - Angle used in picard precheck method
881f40b411bSPeter Brune 
882f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
883f40b411bSPeter Brune 
884f40b411bSPeter Brune    Level: intermediate
885f40b411bSPeter Brune 
8863c7d6663SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetOrder(), SNESLineSearchSetType(), SNESLineSearchSetTolerances(), SNESLineSearchSetDamping(), SNESLineSearchPreCheckPicard()
887f40b411bSPeter Brune @*/
888bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch)
889bf388a1fSBarry Smith {
890bf7f4e0aSPeter Brune   PetscErrorCode    ierr;
8911a4f838cSPeter Brune   const char        *deft = SNESLINESEARCHBASIC;
892bf7f4e0aSPeter Brune   char              type[256];
893bf7f4e0aSPeter Brune   PetscBool         flg, set;
894*dcf2fd19SBarry Smith   PetscViewer       viewer;
895bf388a1fSBarry Smith 
896bf7f4e0aSPeter Brune   PetscFunctionBegin;
8970f51fdf8SToby Isaac   ierr = SNESLineSearchRegisterAll();CHKERRQ(ierr);
898bf7f4e0aSPeter Brune 
899bf7f4e0aSPeter Brune   ierr = PetscObjectOptionsBegin((PetscObject)linesearch);CHKERRQ(ierr);
900f5af7f23SKarl Rupp   if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name;
901a264d7a6SBarry Smith   ierr = PetscOptionsFList("-snes_linesearch_type","Linesearch type","SNESLineSearchSetType",SNESLineSearchList,deft,type,256,&flg);CHKERRQ(ierr);
902bf7f4e0aSPeter Brune   if (flg) {
903f1c6b773SPeter Brune     ierr = SNESLineSearchSetType(linesearch,type);CHKERRQ(ierr);
904bf7f4e0aSPeter Brune   } else if (!((PetscObject)linesearch)->type_name) {
905f1c6b773SPeter Brune     ierr = SNESLineSearchSetType(linesearch,deft);CHKERRQ(ierr);
906bf7f4e0aSPeter Brune   }
907bf7f4e0aSPeter Brune 
908*dcf2fd19SBarry Smith   ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)linesearch),((PetscObject)linesearch)->prefix,"-snes_linesearch_monitor",&viewer,NULL,&set);CHKERRQ(ierr);
909*dcf2fd19SBarry Smith   if (set) {
910*dcf2fd19SBarry Smith     ierr = SNESLineSearchSetDefaultMonitor(linesearch,viewer);CHKERRQ(ierr);
911*dcf2fd19SBarry Smith     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
912*dcf2fd19SBarry Smith   }
913*dcf2fd19SBarry Smith   ierr = SNESLineSearchMonitorSetFromOptions(linesearch,"-snes_linesearch_monitor_solution_update","View correction at each iteration","SNESLineSearchMonitorSolutionUpdate",SNESLineSearchMonitorSolutionUpdate,NULL);CHKERRQ(ierr);
914bf7f4e0aSPeter Brune 
9151a9b3a06SPeter Brune   /* tolerances */
91694ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_minlambda","Minimum step length","SNESLineSearchSetTolerances",linesearch->steptol,&linesearch->steptol,NULL);CHKERRQ(ierr);
91794ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_maxstep","Maximum step size","SNESLineSearchSetTolerances",linesearch->maxstep,&linesearch->maxstep,NULL);CHKERRQ(ierr);
91894ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_rtol","Relative tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->rtol,&linesearch->rtol,NULL);CHKERRQ(ierr);
91994ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_atol","Absolute tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->atol,&linesearch->atol,NULL);CHKERRQ(ierr);
92094ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_ltol","Change in lambda tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->ltol,&linesearch->ltol,NULL);CHKERRQ(ierr);
92194ae4db5SBarry Smith   ierr = PetscOptionsInt("-snes_linesearch_max_it","Maximum iterations for iterative line searches","SNESLineSearchSetTolerances",linesearch->max_its,&linesearch->max_its,NULL);CHKERRQ(ierr);
9227a35526eSPeter Brune 
9231a9b3a06SPeter Brune   /* damping parameters */
92494ae4db5SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_damping","Line search damping and initial step guess","SNESLineSearchSetDamping",linesearch->damping,&linesearch->damping,NULL);CHKERRQ(ierr);
9251a9b3a06SPeter Brune 
92694ae4db5SBarry Smith   ierr = PetscOptionsBool("-snes_linesearch_keeplambda","Use previous lambda as damping","SNESLineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,NULL);CHKERRQ(ierr);
9271a9b3a06SPeter Brune 
9281a9b3a06SPeter Brune   /* precheck */
9297a35526eSPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_precheck_picard","Use a correction that sometimes improves convergence of Picard iteration","SNESLineSearchPreCheckPicard",flg,&flg,&set);CHKERRQ(ierr);
9307a35526eSPeter Brune   if (set) {
9317a35526eSPeter Brune     if (flg) {
9327a35526eSPeter Brune       linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */
933f5af7f23SKarl Rupp 
9347a35526eSPeter Brune       ierr = PetscOptionsReal("-snes_linesearch_precheck_picard_angle","Maximum angle at which to activate the correction",
9350298fd71SBarry Smith                               "none",linesearch->precheck_picard_angle,&linesearch->precheck_picard_angle,NULL);CHKERRQ(ierr);
9367a35526eSPeter Brune       ierr = SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&linesearch->precheck_picard_angle);CHKERRQ(ierr);
9377a35526eSPeter Brune     } else {
9380298fd71SBarry Smith       ierr = SNESLineSearchSetPreCheck(linesearch,NULL,NULL);CHKERRQ(ierr);
9397a35526eSPeter Brune     }
9407a35526eSPeter Brune   }
94194ae4db5SBarry Smith   ierr = PetscOptionsInt("-snes_linesearch_order","Order of approximation used in the line search","SNESLineSearchSetOrder",linesearch->order,&linesearch->order,NULL);CHKERRQ(ierr);
94294ae4db5SBarry Smith   ierr = PetscOptionsBool("-snes_linesearch_norms","Compute final norms in line search","SNESLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,NULL);CHKERRQ(ierr);
9437a35526eSPeter Brune 
9445a9b6599SPeter Brune   if (linesearch->ops->setfromoptions) {
945e55864a3SBarry Smith     (*linesearch->ops->setfromoptions)(PetscOptionsObject,linesearch);CHKERRQ(ierr);
9465a9b6599SPeter Brune   }
9475a9b6599SPeter Brune 
9480633abcbSJed Brown   ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)linesearch);CHKERRQ(ierr);
949bf7f4e0aSPeter Brune   ierr = PetscOptionsEnd();CHKERRQ(ierr);
950bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
951bf7f4e0aSPeter Brune }
952bf7f4e0aSPeter Brune 
953bf7f4e0aSPeter Brune #undef __FUNCT__
954f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchView"
955f40b411bSPeter Brune /*@
956f190f2fcSBarry Smith    SNESLineSearchView - Prints useful information about the line search
957f40b411bSPeter Brune 
958f40b411bSPeter Brune    Input Parameters:
9598e557f58SPeter Brune .  linesearch - linesearch context
960f40b411bSPeter Brune 
961f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
962f40b411bSPeter Brune 
963f40b411bSPeter Brune    Level: intermediate
964f40b411bSPeter Brune 
965*dcf2fd19SBarry Smith .seealso: SNESLineSearchCreate()
966f40b411bSPeter Brune @*/
967bf388a1fSBarry Smith PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer)
968bf388a1fSBarry Smith {
9697f1410a3SPeter Brune   PetscErrorCode ierr;
9707f1410a3SPeter Brune   PetscBool      iascii;
971bf388a1fSBarry Smith 
972bf7f4e0aSPeter Brune   PetscFunctionBegin;
9737f1410a3SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
9747f1410a3SPeter Brune   if (!viewer) {
975ce94432eSBarry Smith     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch),&viewer);CHKERRQ(ierr);
9767f1410a3SPeter Brune   }
9777f1410a3SPeter Brune   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
9787f1410a3SPeter Brune   PetscCheckSameComm(linesearch,1,viewer,2);
979f40b411bSPeter Brune 
980251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
9817f1410a3SPeter Brune   if (iascii) {
982dae58748SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)linesearch,viewer);CHKERRQ(ierr);
9837f1410a3SPeter Brune     if (linesearch->ops->view) {
9847f1410a3SPeter Brune       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
9857f1410a3SPeter Brune       ierr = (*linesearch->ops->view)(linesearch,viewer);CHKERRQ(ierr);
9867f1410a3SPeter Brune       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
9877f1410a3SPeter Brune     }
988c69d1a72SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep,(double)linesearch->steptol);CHKERRQ(ierr);
989c69d1a72SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol,(double)linesearch->atol,(double)linesearch->ltol);CHKERRQ(ierr);
9907f1410a3SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D\n", linesearch->max_its);CHKERRQ(ierr);
9916b2b7091SBarry Smith     if (linesearch->ops->precheck) {
9926b2b7091SBarry Smith       if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) {
9937f1410a3SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  using precheck step to speed up Picard convergence\n", linesearch->max_its);CHKERRQ(ierr);
9947f1410a3SPeter Brune       } else {
9957f1410a3SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  using user-defined precheck step\n", linesearch->max_its);CHKERRQ(ierr);
9967f1410a3SPeter Brune       }
9977f1410a3SPeter Brune     }
9986b2b7091SBarry Smith     if (linesearch->ops->postcheck) {
9997f1410a3SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer,"  using user-defined postcheck step\n", linesearch->max_its);CHKERRQ(ierr);
10007f1410a3SPeter Brune     }
10017f1410a3SPeter Brune   }
1002bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1003bf7f4e0aSPeter Brune }
1004bf7f4e0aSPeter Brune 
1005bf7f4e0aSPeter Brune #undef __FUNCT__
1006f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetType"
1007ea5d4fccSPeter Brune /*@C
1008f1c6b773SPeter Brune    SNESLineSearchSetType - Sets the linesearch type
1009f40b411bSPeter Brune 
1010f190f2fcSBarry Smith    Logically Collective on SNESLineSearch
1011f190f2fcSBarry Smith 
1012f40b411bSPeter Brune    Input Parameters:
10138e557f58SPeter Brune +  linesearch - linesearch context
1014f40b411bSPeter Brune -  type - The type of line search to be used
1015f40b411bSPeter Brune 
1016cd7522eaSPeter Brune    Available Types:
1017cd7522eaSPeter Brune +  basic - Simple damping line search.
10188e557f58SPeter Brune .  bt - Backtracking line search over the L2 norm of the function
10198e557f58SPeter Brune .  l2 - Secant line search over the L2 norm of the function
10208e557f58SPeter Brune .  cp - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x)
1021d4c6564cSPatrick Farrell .  nleqerr - Affine-covariant error-oriented linesearch
10228e557f58SPeter Brune -  shell - User provided SNESLineSearch implementation
1023cd7522eaSPeter Brune 
1024f40b411bSPeter Brune    Level: intermediate
1025f40b411bSPeter Brune 
1026f1c6b773SPeter Brune .seealso: SNESLineSearchCreate()
1027f40b411bSPeter Brune @*/
102819fd82e9SBarry Smith PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type)
1029bf7f4e0aSPeter Brune {
1030f1c6b773SPeter Brune   PetscErrorCode ierr,(*r)(SNESLineSearch);
1031bf7f4e0aSPeter Brune   PetscBool      match;
1032bf7f4e0aSPeter Brune 
1033bf7f4e0aSPeter Brune   PetscFunctionBegin;
1034f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1035bf7f4e0aSPeter Brune   PetscValidCharPointer(type,2);
1036bf7f4e0aSPeter Brune 
1037251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)linesearch,type,&match);CHKERRQ(ierr);
1038bf7f4e0aSPeter Brune   if (match) PetscFunctionReturn(0);
1039bf7f4e0aSPeter Brune 
10401c9cd337SJed Brown   ierr = PetscFunctionListFind(SNESLineSearchList,type,&r);CHKERRQ(ierr);
1041bf7f4e0aSPeter Brune   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type);
1042bf7f4e0aSPeter Brune   /* Destroy the previous private linesearch context */
1043bf7f4e0aSPeter Brune   if (linesearch->ops->destroy) {
1044bf7f4e0aSPeter Brune     ierr = (*(linesearch)->ops->destroy)(linesearch);CHKERRQ(ierr);
1045f5af7f23SKarl Rupp 
10460298fd71SBarry Smith     linesearch->ops->destroy = NULL;
1047bf7f4e0aSPeter Brune   }
1048f1c6b773SPeter Brune   /* Reinitialize function pointers in SNESLineSearchOps structure */
1049bf7f4e0aSPeter Brune   linesearch->ops->apply          = 0;
1050bf7f4e0aSPeter Brune   linesearch->ops->view           = 0;
1051bf7f4e0aSPeter Brune   linesearch->ops->setfromoptions = 0;
1052bf7f4e0aSPeter Brune   linesearch->ops->destroy        = 0;
1053bf7f4e0aSPeter Brune 
1054bf7f4e0aSPeter Brune   ierr = PetscObjectChangeTypeName((PetscObject)linesearch,type);CHKERRQ(ierr);
1055bf7f4e0aSPeter Brune   ierr = (*r)(linesearch);CHKERRQ(ierr);
1056bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1057bf7f4e0aSPeter Brune }
1058bf7f4e0aSPeter Brune 
1059bf7f4e0aSPeter Brune #undef __FUNCT__
1060f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetSNES"
1061f40b411bSPeter Brune /*@
106278bcb3b5SPeter Brune    SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation.
1063f40b411bSPeter Brune 
1064f40b411bSPeter Brune    Input Parameters:
10658e557f58SPeter Brune +  linesearch - linesearch context
1066f40b411bSPeter Brune -  snes - The snes instance
1067f40b411bSPeter Brune 
106878bcb3b5SPeter Brune    Level: developer
106978bcb3b5SPeter Brune 
107078bcb3b5SPeter Brune    Notes:
1071f190f2fcSBarry Smith    This happens automatically when the line search is obtained/created with
10727601faf0SJed Brown    SNESGetLineSearch().  This routine is therefore mainly called within SNES
107378bcb3b5SPeter Brune    implementations.
1074f40b411bSPeter Brune 
10758141a3b9SPeter Brune    Level: developer
1076f40b411bSPeter Brune 
1077cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
1078f40b411bSPeter Brune @*/
1079bf388a1fSBarry Smith PetscErrorCode  SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes)
1080bf388a1fSBarry Smith {
1081bf7f4e0aSPeter Brune   PetscFunctionBegin;
1082f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1083bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,2);
1084bf7f4e0aSPeter Brune   linesearch->snes = snes;
1085bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1086bf7f4e0aSPeter Brune }
1087bf7f4e0aSPeter Brune 
1088bf7f4e0aSPeter Brune #undef __FUNCT__
1089f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetSNES"
1090f40b411bSPeter Brune /*@
10918141a3b9SPeter Brune    SNESLineSearchGetSNES - Gets the SNES instance associated with the line search.
10928141a3b9SPeter Brune    Having an associated SNES is necessary because most line search implementations must be able to
10938141a3b9SPeter Brune    evaluate the function using SNESComputeFunction() for the associated SNES.  This routine
10948141a3b9SPeter Brune    is used in the line search implementations when one must get this associated SNES instance.
1095f40b411bSPeter Brune 
1096f40b411bSPeter Brune    Input Parameters:
10978e557f58SPeter Brune .  linesearch - linesearch context
1098f40b411bSPeter Brune 
1099f40b411bSPeter Brune    Output Parameters:
1100f40b411bSPeter Brune .  snes - The snes instance
1101f40b411bSPeter Brune 
11028141a3b9SPeter Brune    Level: developer
1103f40b411bSPeter Brune 
1104cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
1105f40b411bSPeter Brune @*/
1106bf388a1fSBarry Smith PetscErrorCode  SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes)
1107bf388a1fSBarry Smith {
1108bf7f4e0aSPeter Brune   PetscFunctionBegin;
1109f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
11106a388c36SPeter Brune   PetscValidPointer(snes, 2);
1111bf7f4e0aSPeter Brune   *snes = linesearch->snes;
1112bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1113bf7f4e0aSPeter Brune }
1114bf7f4e0aSPeter Brune 
11156a388c36SPeter Brune #undef __FUNCT__
1116f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetLambda"
1117f40b411bSPeter Brune /*@
1118f1c6b773SPeter Brune    SNESLineSearchGetLambda - Gets the last linesearch steplength discovered.
1119f40b411bSPeter Brune 
1120f40b411bSPeter Brune    Input Parameters:
11218e557f58SPeter Brune .  linesearch - linesearch context
1122f40b411bSPeter Brune 
1123f40b411bSPeter Brune    Output Parameters:
1124cd7522eaSPeter Brune .  lambda - The last steplength computed during SNESLineSearchApply()
1125f40b411bSPeter Brune 
112678bcb3b5SPeter Brune    Level: advanced
112778bcb3b5SPeter Brune 
11288e557f58SPeter Brune    Notes:
11298e557f58SPeter Brune    This is useful in methods where the solver is ill-scaled and
113078bcb3b5SPeter Brune    requires some adaptive notion of the difference in scale between the
113178bcb3b5SPeter Brune    solution and the function.  For instance, SNESQN may be scaled by the
113278bcb3b5SPeter Brune    line search lambda using the argument -snes_qn_scaling ls.
113378bcb3b5SPeter Brune 
1134cd7522eaSPeter Brune .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetDamping(), SNESLineSearchApply()
1135f40b411bSPeter Brune @*/
1136f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda)
11376a388c36SPeter Brune {
11386a388c36SPeter Brune   PetscFunctionBegin;
1139f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
11406a388c36SPeter Brune   PetscValidPointer(lambda, 2);
11416a388c36SPeter Brune   *lambda = linesearch->lambda;
11426a388c36SPeter Brune   PetscFunctionReturn(0);
11436a388c36SPeter Brune }
11446a388c36SPeter Brune 
11456a388c36SPeter Brune #undef __FUNCT__
1146f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetLambda"
1147f40b411bSPeter Brune /*@
1148f1c6b773SPeter Brune    SNESLineSearchSetLambda - Sets the linesearch steplength.
1149f40b411bSPeter Brune 
1150f40b411bSPeter Brune    Input Parameters:
11518e557f58SPeter Brune +  linesearch - linesearch context
1152f40b411bSPeter Brune -  lambda - The last steplength.
1153f40b411bSPeter Brune 
1154cd7522eaSPeter Brune    Notes:
1155f190f2fcSBarry Smith    This routine is typically used within implementations of SNESLineSearchApply()
1156cd7522eaSPeter Brune    to set the final steplength.  This routine (and SNESLineSearchGetLambda()) were
1157cd7522eaSPeter Brune    added in order to facilitate Quasi-Newton methods that use the previous steplength
1158cd7522eaSPeter Brune    as an inner scaling parameter.
1159cd7522eaSPeter Brune 
116078bcb3b5SPeter Brune    Level: advanced
1161f40b411bSPeter Brune 
1162f1c6b773SPeter Brune .seealso: SNESLineSearchGetLambda()
1163f40b411bSPeter Brune @*/
1164f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
11656a388c36SPeter Brune {
11666a388c36SPeter Brune   PetscFunctionBegin;
1167f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
11686a388c36SPeter Brune   linesearch->lambda = lambda;
11696a388c36SPeter Brune   PetscFunctionReturn(0);
11706a388c36SPeter Brune }
11716a388c36SPeter Brune 
11726a388c36SPeter Brune #undef  __FUNCT__
1173f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetTolerances"
1174f40b411bSPeter Brune /*@
11753c7d6663SPeter Brune    SNESLineSearchGetTolerances - Gets the tolerances for the linesearch.  These include
117678bcb3b5SPeter Brune    tolerances for the relative and absolute change in the function norm, the change
117778bcb3b5SPeter Brune    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
117878bcb3b5SPeter Brune    and the maximum number of iterations the line search procedure may take.
1179f40b411bSPeter Brune 
1180f40b411bSPeter Brune    Input Parameters:
11818e557f58SPeter Brune .  linesearch - linesearch context
1182f40b411bSPeter Brune 
1183f40b411bSPeter Brune    Output Parameters:
1184516fe3c3SPeter Brune +  steptol - The minimum steplength
11856cc8e53bSPeter Brune .  maxstep - The maximum steplength
1186516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
1187516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
1188516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
1189516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
1190f40b411bSPeter Brune 
119178bcb3b5SPeter Brune    Level: intermediate
119278bcb3b5SPeter Brune 
119378bcb3b5SPeter Brune    Notes:
119478bcb3b5SPeter Brune    Different line searches may implement these parameters slightly differently as
11953c7d6663SPeter Brune    the type requires.
1196516fe3c3SPeter Brune 
1197f1c6b773SPeter Brune .seealso: SNESLineSearchSetTolerances()
1198f40b411bSPeter Brune @*/
1199f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
12006a388c36SPeter Brune {
12016a388c36SPeter Brune   PetscFunctionBegin;
1202f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1203516fe3c3SPeter Brune   if (steptol) {
12046a388c36SPeter Brune     PetscValidPointer(steptol, 2);
12056a388c36SPeter Brune     *steptol = linesearch->steptol;
1206516fe3c3SPeter Brune   }
1207516fe3c3SPeter Brune   if (maxstep) {
1208516fe3c3SPeter Brune     PetscValidPointer(maxstep, 3);
1209516fe3c3SPeter Brune     *maxstep = linesearch->maxstep;
1210516fe3c3SPeter Brune   }
1211516fe3c3SPeter Brune   if (rtol) {
1212516fe3c3SPeter Brune     PetscValidPointer(rtol, 4);
1213516fe3c3SPeter Brune     *rtol = linesearch->rtol;
1214516fe3c3SPeter Brune   }
1215516fe3c3SPeter Brune   if (atol) {
1216516fe3c3SPeter Brune     PetscValidPointer(atol, 5);
1217516fe3c3SPeter Brune     *atol = linesearch->atol;
1218516fe3c3SPeter Brune   }
1219516fe3c3SPeter Brune   if (ltol) {
1220516fe3c3SPeter Brune     PetscValidPointer(ltol, 6);
1221516fe3c3SPeter Brune     *ltol = linesearch->ltol;
1222516fe3c3SPeter Brune   }
1223516fe3c3SPeter Brune   if (max_its) {
1224516fe3c3SPeter Brune     PetscValidPointer(max_its, 7);
1225516fe3c3SPeter Brune     *max_its = linesearch->max_its;
1226516fe3c3SPeter Brune   }
12276a388c36SPeter Brune   PetscFunctionReturn(0);
12286a388c36SPeter Brune }
12296a388c36SPeter Brune 
12306a388c36SPeter Brune #undef  __FUNCT__
1231f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetTolerances"
1232f40b411bSPeter Brune /*@
12333c7d6663SPeter Brune    SNESLineSearchSetTolerances -  Gets the tolerances for the linesearch.  These include
123478bcb3b5SPeter Brune    tolerances for the relative and absolute change in the function norm, the change
123578bcb3b5SPeter Brune    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
123678bcb3b5SPeter Brune    and the maximum number of iterations the line search procedure may take.
1237f40b411bSPeter Brune 
1238f40b411bSPeter Brune    Input Parameters:
12398e557f58SPeter Brune +  linesearch - linesearch context
1240516fe3c3SPeter Brune .  steptol - The minimum steplength
12416cc8e53bSPeter Brune .  maxstep - The maximum steplength
1242516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
1243516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
1244516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
1245516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
1246f40b411bSPeter Brune 
124778bcb3b5SPeter Brune    Notes:
12483c7d6663SPeter Brune    The user may choose to not set any of the tolerances using PETSC_DEFAULT in
124978bcb3b5SPeter Brune    place of an argument.
1250f40b411bSPeter Brune 
125178bcb3b5SPeter Brune    Level: intermediate
1252516fe3c3SPeter Brune 
1253f1c6b773SPeter Brune .seealso: SNESLineSearchGetTolerances()
1254f40b411bSPeter Brune @*/
1255f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its)
12566a388c36SPeter Brune {
12576a388c36SPeter Brune   PetscFunctionBegin;
1258f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1259d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,steptol,2);
1260d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,maxstep,3);
1261d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,rtol,4);
1262d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,atol,5);
1263d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,ltol,6);
1264d3952184SSatish Balay   PetscValidLogicalCollectiveInt(linesearch,max_its,7);
1265d3952184SSatish Balay 
1266d3952184SSatish Balay   if (steptol!= PETSC_DEFAULT) {
1267ce94432eSBarry Smith     if (steptol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Minimum step length %14.12e must be non-negative",(double)steptol);
12686a388c36SPeter Brune     linesearch->steptol = steptol;
1269d3952184SSatish Balay   }
1270d3952184SSatish Balay 
1271d3952184SSatish Balay   if (maxstep!= PETSC_DEFAULT) {
1272ce94432eSBarry Smith     if (maxstep < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum step length %14.12e must be non-negative",(double)maxstep);
1273516fe3c3SPeter Brune     linesearch->maxstep = maxstep;
1274d3952184SSatish Balay   }
1275d3952184SSatish Balay 
1276d3952184SSatish Balay   if (rtol != PETSC_DEFAULT) {
1277ce94432eSBarry 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);
1278516fe3c3SPeter Brune     linesearch->rtol = rtol;
1279d3952184SSatish Balay   }
1280d3952184SSatish Balay 
1281d3952184SSatish Balay   if (atol != PETSC_DEFAULT) {
1282ce94432eSBarry Smith     if (atol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %14.12e must be non-negative",(double)atol);
1283516fe3c3SPeter Brune     linesearch->atol = atol;
1284d3952184SSatish Balay   }
1285d3952184SSatish Balay 
1286d3952184SSatish Balay   if (ltol != PETSC_DEFAULT) {
1287ce94432eSBarry Smith     if (ltol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Labmda tolerance %14.12e must be non-negative",(double)ltol);
1288516fe3c3SPeter Brune     linesearch->ltol = ltol;
1289d3952184SSatish Balay   }
1290d3952184SSatish Balay 
1291d3952184SSatish Balay   if (max_its != PETSC_DEFAULT) {
1292ce94432eSBarry Smith     if (max_its < 0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",max_its);
1293516fe3c3SPeter Brune     linesearch->max_its = max_its;
1294d3952184SSatish Balay   }
12956a388c36SPeter Brune   PetscFunctionReturn(0);
12966a388c36SPeter Brune }
12976a388c36SPeter Brune 
12986a388c36SPeter Brune #undef __FUNCT__
1299f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetDamping"
1300f40b411bSPeter Brune /*@
1301f1c6b773SPeter Brune    SNESLineSearchGetDamping - Gets the line search damping parameter.
1302f40b411bSPeter Brune 
1303f40b411bSPeter Brune    Input Parameters:
13048e557f58SPeter Brune .  linesearch - linesearch context
1305f40b411bSPeter Brune 
1306f40b411bSPeter Brune    Output Parameters:
13078e557f58SPeter Brune .  damping - The damping parameter
1308f40b411bSPeter Brune 
130978bcb3b5SPeter Brune    Level: advanced
1310f40b411bSPeter Brune 
131178bcb3b5SPeter Brune .seealso: SNESLineSearchGetStepTolerance(), SNESQN
1312f40b411bSPeter Brune @*/
1313f40b411bSPeter Brune 
1314f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping)
13156a388c36SPeter Brune {
13166a388c36SPeter Brune   PetscFunctionBegin;
1317f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
13186a388c36SPeter Brune   PetscValidPointer(damping, 2);
13196a388c36SPeter Brune   *damping = linesearch->damping;
13206a388c36SPeter Brune   PetscFunctionReturn(0);
13216a388c36SPeter Brune }
13226a388c36SPeter Brune 
13236a388c36SPeter Brune #undef __FUNCT__
1324f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetDamping"
1325f40b411bSPeter Brune /*@
1326f1c6b773SPeter Brune    SNESLineSearchSetDamping - Sets the line search damping paramter.
1327f40b411bSPeter Brune 
1328f40b411bSPeter Brune    Input Parameters:
132903fd4120SBarry Smith +  linesearch - linesearch context
133003fd4120SBarry Smith -  damping - The damping parameter
1331f40b411bSPeter Brune 
133203fd4120SBarry Smith    Options Database:
133303fd4120SBarry Smith .   -snes_linesearch_damping
1334f40b411bSPeter Brune    Level: intermediate
1335f40b411bSPeter Brune 
1336cd7522eaSPeter Brune    Notes:
1337cd7522eaSPeter Brune    The basic line search merely takes the update step scaled by the damping parameter.
1338cd7522eaSPeter Brune    The use of the damping parameter in the l2 and cp line searches is much more subtle;
133978bcb3b5SPeter Brune    it is used as a starting point in calculating the secant step. However, the eventual
1340cd7522eaSPeter Brune    step may be of greater length than the damping parameter.  In the bt line search it is
1341cd7522eaSPeter Brune    used as the maximum possible step length, as the bt line search only backtracks.
1342cd7522eaSPeter Brune 
1343f1c6b773SPeter Brune .seealso: SNESLineSearchGetDamping()
1344f40b411bSPeter Brune @*/
1345f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping)
13466a388c36SPeter Brune {
13476a388c36SPeter Brune   PetscFunctionBegin;
1348f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
13496a388c36SPeter Brune   linesearch->damping = damping;
13506a388c36SPeter Brune   PetscFunctionReturn(0);
13516a388c36SPeter Brune }
13526a388c36SPeter Brune 
13536a388c36SPeter Brune #undef __FUNCT__
135459405d5eSPeter Brune #define __FUNCT__ "SNESLineSearchGetOrder"
135559405d5eSPeter Brune /*@
135659405d5eSPeter Brune    SNESLineSearchGetOrder - Gets the line search approximation order.
135759405d5eSPeter Brune 
135859405d5eSPeter Brune    Input Parameters:
135978bcb3b5SPeter Brune .  linesearch - linesearch context
136059405d5eSPeter Brune 
136159405d5eSPeter Brune    Output Parameters:
13628e557f58SPeter Brune .  order - The order
136359405d5eSPeter Brune 
136478bcb3b5SPeter Brune    Possible Values for order:
13653c7d6663SPeter Brune +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
13663c7d6663SPeter Brune .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
13673c7d6663SPeter Brune -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
136878bcb3b5SPeter Brune 
136959405d5eSPeter Brune    Level: intermediate
137059405d5eSPeter Brune 
137159405d5eSPeter Brune .seealso: SNESLineSearchSetOrder()
137259405d5eSPeter Brune @*/
137359405d5eSPeter Brune 
1374b000cd8dSPeter Brune PetscErrorCode  SNESLineSearchGetOrder(SNESLineSearch linesearch,PetscInt *order)
137559405d5eSPeter Brune {
137659405d5eSPeter Brune   PetscFunctionBegin;
137759405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
137859405d5eSPeter Brune   PetscValidPointer(order, 2);
137959405d5eSPeter Brune   *order = linesearch->order;
138059405d5eSPeter Brune   PetscFunctionReturn(0);
138159405d5eSPeter Brune }
138259405d5eSPeter Brune 
138359405d5eSPeter Brune #undef __FUNCT__
138459405d5eSPeter Brune #define __FUNCT__ "SNESLineSearchSetOrder"
138559405d5eSPeter Brune /*@
138659405d5eSPeter Brune    SNESLineSearchSetOrder - Sets the line search damping paramter.
138759405d5eSPeter Brune 
138859405d5eSPeter Brune    Input Parameters:
138978bcb3b5SPeter Brune .  linesearch - linesearch context
139078bcb3b5SPeter Brune .  order - The damping parameter
139159405d5eSPeter Brune 
139259405d5eSPeter Brune    Level: intermediate
139359405d5eSPeter Brune 
139478bcb3b5SPeter Brune    Possible Values for order:
13953c7d6663SPeter Brune +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
13963c7d6663SPeter Brune .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
13973c7d6663SPeter Brune -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
139878bcb3b5SPeter Brune 
139959405d5eSPeter Brune    Notes:
140059405d5eSPeter Brune    Variable orders are supported by the following line searches:
140178bcb3b5SPeter Brune +  bt - cubic and quadratic
140278bcb3b5SPeter Brune -  cp - linear and quadratic
140359405d5eSPeter Brune 
140459405d5eSPeter Brune .seealso: SNESLineSearchGetOrder()
140559405d5eSPeter Brune @*/
1406b000cd8dSPeter Brune PetscErrorCode  SNESLineSearchSetOrder(SNESLineSearch linesearch,PetscInt order)
140759405d5eSPeter Brune {
140859405d5eSPeter Brune   PetscFunctionBegin;
140959405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
141059405d5eSPeter Brune   linesearch->order = order;
141159405d5eSPeter Brune   PetscFunctionReturn(0);
141259405d5eSPeter Brune }
141359405d5eSPeter Brune 
141459405d5eSPeter Brune #undef __FUNCT__
1415f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetNorms"
1416f40b411bSPeter Brune /*@
1417f1c6b773SPeter Brune    SNESLineSearchGetNorms - Gets the norms for for X, Y, and F.
1418f40b411bSPeter Brune 
1419f40b411bSPeter Brune    Input Parameters:
142078bcb3b5SPeter Brune .  linesearch - linesearch context
1421f40b411bSPeter Brune 
1422f40b411bSPeter Brune    Output Parameters:
1423f40b411bSPeter Brune +  xnorm - The norm of the current solution
1424f40b411bSPeter Brune .  fnorm - The norm of the current function
1425f40b411bSPeter Brune -  ynorm - The norm of the current update
1426f40b411bSPeter Brune 
1427cd7522eaSPeter Brune    Notes:
1428cd7522eaSPeter Brune    This function is mainly called from SNES implementations.
1429cd7522eaSPeter Brune 
143078bcb3b5SPeter Brune    Level: developer
1431f40b411bSPeter Brune 
1432f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms() SNESLineSearchGetVecs()
1433f40b411bSPeter Brune @*/
1434f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm)
1435bf7f4e0aSPeter Brune {
1436bf7f4e0aSPeter Brune   PetscFunctionBegin;
1437f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1438f5af7f23SKarl Rupp   if (xnorm) *xnorm = linesearch->xnorm;
1439f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
1440f5af7f23SKarl Rupp   if (ynorm) *ynorm = linesearch->ynorm;
1441bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1442bf7f4e0aSPeter Brune }
1443bf7f4e0aSPeter Brune 
1444e7058c64SPeter Brune #undef __FUNCT__
1445f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetNorms"
1446f40b411bSPeter Brune /*@
1447f1c6b773SPeter Brune    SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F.
1448f40b411bSPeter Brune 
1449f40b411bSPeter Brune    Input Parameters:
145078bcb3b5SPeter Brune +  linesearch - linesearch context
1451f40b411bSPeter Brune .  xnorm - The norm of the current solution
1452f40b411bSPeter Brune .  fnorm - The norm of the current function
1453f40b411bSPeter Brune -  ynorm - The norm of the current update
1454f40b411bSPeter Brune 
145578bcb3b5SPeter Brune    Level: advanced
1456f40b411bSPeter Brune 
1457f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1458f40b411bSPeter Brune @*/
1459f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
14606a388c36SPeter Brune {
14616a388c36SPeter Brune   PetscFunctionBegin;
1462f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
14636a388c36SPeter Brune   linesearch->xnorm = xnorm;
14646a388c36SPeter Brune   linesearch->fnorm = fnorm;
14656a388c36SPeter Brune   linesearch->ynorm = ynorm;
14666a388c36SPeter Brune   PetscFunctionReturn(0);
14676a388c36SPeter Brune }
14686a388c36SPeter Brune 
14696a388c36SPeter Brune #undef __FUNCT__
1470f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchComputeNorms"
1471f40b411bSPeter Brune /*@
1472f1c6b773SPeter Brune    SNESLineSearchComputeNorms - Computes the norms of X, F, and Y.
1473f40b411bSPeter Brune 
1474f40b411bSPeter Brune    Input Parameters:
147578bcb3b5SPeter Brune .  linesearch - linesearch context
1476f40b411bSPeter Brune 
1477f40b411bSPeter Brune    Options Database Keys:
14788e557f58SPeter Brune .   -snes_linesearch_norms - turn norm computation on or off
1479f40b411bSPeter Brune 
1480f40b411bSPeter Brune    Level: intermediate
1481f40b411bSPeter Brune 
148278bcb3b5SPeter Brune .seealso: SNESLineSearchGetNorms, SNESLineSearchSetNorms(), SNESLineSearchSetComputeNorms()
1483f40b411bSPeter Brune @*/
1484f1c6b773SPeter Brune PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
14856a388c36SPeter Brune {
14866a388c36SPeter Brune   PetscErrorCode ierr;
14879bd66eb0SPeter Brune   SNES           snes;
1488bf388a1fSBarry Smith 
14896a388c36SPeter Brune   PetscFunctionBegin;
14906a388c36SPeter Brune   if (linesearch->norms) {
14919bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
1492f1c6b773SPeter Brune       ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
14939bd66eb0SPeter Brune       ierr = VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
14949bd66eb0SPeter Brune       ierr = VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
14959bd66eb0SPeter Brune       ierr = (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);CHKERRQ(ierr);
14969bd66eb0SPeter Brune     } else {
14976a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_func,   NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
14986a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_sol,    NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
14996a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
15006a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_func,     NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
15016a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_sol,      NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
15026a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_update,   NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
15036a388c36SPeter Brune     }
15049bd66eb0SPeter Brune   }
15056a388c36SPeter Brune   PetscFunctionReturn(0);
15066a388c36SPeter Brune }
15076a388c36SPeter Brune 
15086f263ca3SPeter Brune #undef __FUNCT__
15096f263ca3SPeter Brune #define __FUNCT__ "SNESLineSearchSetComputeNorms"
15106f263ca3SPeter Brune /*@
15116f263ca3SPeter Brune    SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search.
15126f263ca3SPeter Brune 
15136f263ca3SPeter Brune    Input Parameters:
151478bcb3b5SPeter Brune +  linesearch  - linesearch context
151578bcb3b5SPeter Brune -  flg  - indicates whether or not to compute norms
15166f263ca3SPeter Brune 
15176f263ca3SPeter Brune    Options Database Keys:
15188e557f58SPeter Brune .   -snes_linesearch_norms - turn norm computation on or off
15196f263ca3SPeter Brune 
15206f263ca3SPeter Brune    Notes:
15211a4f838cSPeter Brune    This is most relevant to the SNESLINESEARCHBASIC line search type.
15226f263ca3SPeter Brune 
15236f263ca3SPeter Brune    Level: intermediate
15246f263ca3SPeter Brune 
15251a4f838cSPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetNorms(), SNESLineSearchComputeNorms(), SNESLINESEARCHBASIC
15266f263ca3SPeter Brune @*/
15276f263ca3SPeter Brune PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg)
15286f263ca3SPeter Brune {
15296f263ca3SPeter Brune   PetscFunctionBegin;
15306f263ca3SPeter Brune   linesearch->norms = flg;
15316f263ca3SPeter Brune   PetscFunctionReturn(0);
15326f263ca3SPeter Brune }
15336f263ca3SPeter Brune 
15346a388c36SPeter Brune #undef __FUNCT__
1535f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetVecs"
1536f40b411bSPeter Brune /*@
1537f1c6b773SPeter Brune    SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context
1538f40b411bSPeter Brune 
1539f40b411bSPeter Brune    Input Parameters:
154078bcb3b5SPeter Brune .  linesearch - linesearch context
1541f40b411bSPeter Brune 
1542f40b411bSPeter Brune    Output Parameters:
15436232e825SPeter Brune +  X - Solution vector
15446232e825SPeter Brune .  F - Function vector
15456232e825SPeter Brune .  Y - Search direction vector
15466232e825SPeter Brune .  W - Solution work vector
15476232e825SPeter Brune -  G - Function work vector
15486232e825SPeter Brune 
15497bba9028SPeter Brune    Notes:
15507bba9028SPeter Brune    At the beginning of a line search application, X should contain a
15516232e825SPeter Brune    solution and the vector F the function computed at X.  At the end of the
15526232e825SPeter Brune    line search application, X should contain the new solution, and F the
15536232e825SPeter Brune    function evaluated at the new solution.
1554f40b411bSPeter Brune 
15552a7a6963SBarry Smith    These vectors are owned by the SNESLineSearch and should not be destroyed by the caller
15562a7a6963SBarry Smith 
155778bcb3b5SPeter Brune    Level: advanced
1558f40b411bSPeter Brune 
1559f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1560f40b411bSPeter Brune @*/
1561bf388a1fSBarry Smith PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G)
1562bf388a1fSBarry Smith {
15636a388c36SPeter Brune   PetscFunctionBegin;
1564f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
15656a388c36SPeter Brune   if (X) {
15666a388c36SPeter Brune     PetscValidPointer(X, 2);
15676a388c36SPeter Brune     *X = linesearch->vec_sol;
15686a388c36SPeter Brune   }
15696a388c36SPeter Brune   if (F) {
15706a388c36SPeter Brune     PetscValidPointer(F, 3);
15716a388c36SPeter Brune     *F = linesearch->vec_func;
15726a388c36SPeter Brune   }
15736a388c36SPeter Brune   if (Y) {
15746a388c36SPeter Brune     PetscValidPointer(Y, 4);
15756a388c36SPeter Brune     *Y = linesearch->vec_update;
15766a388c36SPeter Brune   }
15776a388c36SPeter Brune   if (W) {
15786a388c36SPeter Brune     PetscValidPointer(W, 5);
15796a388c36SPeter Brune     *W = linesearch->vec_sol_new;
15806a388c36SPeter Brune   }
15816a388c36SPeter Brune   if (G) {
15826a388c36SPeter Brune     PetscValidPointer(G, 6);
15836a388c36SPeter Brune     *G = linesearch->vec_func_new;
15846a388c36SPeter Brune   }
15856a388c36SPeter Brune   PetscFunctionReturn(0);
15866a388c36SPeter Brune }
15876a388c36SPeter Brune 
15886a388c36SPeter Brune #undef __FUNCT__
1589f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetVecs"
1590f40b411bSPeter Brune /*@
1591f1c6b773SPeter Brune    SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context
1592f40b411bSPeter Brune 
1593f40b411bSPeter Brune    Input Parameters:
159478bcb3b5SPeter Brune +  linesearch - linesearch context
15956232e825SPeter Brune .  X - Solution vector
15966232e825SPeter Brune .  F - Function vector
15976232e825SPeter Brune .  Y - Search direction vector
15986232e825SPeter Brune .  W - Solution work vector
15996232e825SPeter Brune -  G - Function work vector
1600f40b411bSPeter Brune 
160178bcb3b5SPeter Brune    Level: advanced
1602f40b411bSPeter Brune 
1603f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms(), SNESLineSearchGetVecs()
1604f40b411bSPeter Brune @*/
1605bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G)
1606bf388a1fSBarry Smith {
16076a388c36SPeter Brune   PetscFunctionBegin;
1608f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
16096a388c36SPeter Brune   if (X) {
16106a388c36SPeter Brune     PetscValidHeaderSpecific(X,VEC_CLASSID,2);
16116a388c36SPeter Brune     linesearch->vec_sol = X;
16126a388c36SPeter Brune   }
16136a388c36SPeter Brune   if (F) {
16146a388c36SPeter Brune     PetscValidHeaderSpecific(F,VEC_CLASSID,3);
16156a388c36SPeter Brune     linesearch->vec_func = F;
16166a388c36SPeter Brune   }
16176a388c36SPeter Brune   if (Y) {
16186a388c36SPeter Brune     PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
16196a388c36SPeter Brune     linesearch->vec_update = Y;
16206a388c36SPeter Brune   }
16216a388c36SPeter Brune   if (W) {
16226a388c36SPeter Brune     PetscValidHeaderSpecific(W,VEC_CLASSID,5);
16236a388c36SPeter Brune     linesearch->vec_sol_new = W;
16246a388c36SPeter Brune   }
16256a388c36SPeter Brune   if (G) {
16266a388c36SPeter Brune     PetscValidHeaderSpecific(G,VEC_CLASSID,6);
16276a388c36SPeter Brune     linesearch->vec_func_new = G;
16286a388c36SPeter Brune   }
16296a388c36SPeter Brune   PetscFunctionReturn(0);
16306a388c36SPeter Brune }
16316a388c36SPeter Brune 
16326a388c36SPeter Brune #undef __FUNCT__
1633f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchAppendOptionsPrefix"
1634e7058c64SPeter Brune /*@C
1635f1c6b773SPeter Brune    SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1636e7058c64SPeter Brune    SNES options in the database.
1637e7058c64SPeter Brune 
1638cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
1639e7058c64SPeter Brune 
1640e7058c64SPeter Brune    Input Parameters:
1641e7058c64SPeter Brune +  snes - the SNES context
1642e7058c64SPeter Brune -  prefix - the prefix to prepend to all option names
1643e7058c64SPeter Brune 
1644e7058c64SPeter Brune    Notes:
1645e7058c64SPeter Brune    A hyphen (-) must NOT be given at the beginning of the prefix name.
1646e7058c64SPeter Brune    The first character of all runtime options is AUTOMATICALLY the hyphen.
1647e7058c64SPeter Brune 
1648e7058c64SPeter Brune    Level: advanced
1649e7058c64SPeter Brune 
1650f1c6b773SPeter Brune .keywords: SNESLineSearch, append, options, prefix, database
1651e7058c64SPeter Brune 
1652e7058c64SPeter Brune .seealso: SNESGetOptionsPrefix()
1653e7058c64SPeter Brune @*/
1654f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[])
1655e7058c64SPeter Brune {
1656e7058c64SPeter Brune   PetscErrorCode ierr;
1657e7058c64SPeter Brune 
1658e7058c64SPeter Brune   PetscFunctionBegin;
1659f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1660e7058c64SPeter Brune   ierr = PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1661e7058c64SPeter Brune   PetscFunctionReturn(0);
1662e7058c64SPeter Brune }
1663e7058c64SPeter Brune 
1664e7058c64SPeter Brune #undef __FUNCT__
1665f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetOptionsPrefix"
1666e7058c64SPeter Brune /*@C
1667f1c6b773SPeter Brune    SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all
1668f1c6b773SPeter Brune    SNESLineSearch options in the database.
1669e7058c64SPeter Brune 
1670e7058c64SPeter Brune    Not Collective
1671e7058c64SPeter Brune 
1672e7058c64SPeter Brune    Input Parameter:
1673cd7522eaSPeter Brune .  linesearch - the SNESLineSearch context
1674e7058c64SPeter Brune 
1675e7058c64SPeter Brune    Output Parameter:
1676e7058c64SPeter Brune .  prefix - pointer to the prefix string used
1677e7058c64SPeter Brune 
16788e557f58SPeter Brune    Notes:
16798e557f58SPeter Brune    On the fortran side, the user should pass in a string 'prefix' of
1680e7058c64SPeter Brune    sufficient length to hold the prefix.
1681e7058c64SPeter Brune 
1682e7058c64SPeter Brune    Level: advanced
1683e7058c64SPeter Brune 
1684f1c6b773SPeter Brune .keywords: SNESLineSearch, get, options, prefix, database
1685e7058c64SPeter Brune 
1686e7058c64SPeter Brune .seealso: SNESAppendOptionsPrefix()
1687e7058c64SPeter Brune @*/
1688f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[])
1689e7058c64SPeter Brune {
1690e7058c64SPeter Brune   PetscErrorCode ierr;
1691e7058c64SPeter Brune 
1692e7058c64SPeter Brune   PetscFunctionBegin;
1693f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1694e7058c64SPeter Brune   ierr = PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1695e7058c64SPeter Brune   PetscFunctionReturn(0);
1696e7058c64SPeter Brune }
1697bf7f4e0aSPeter Brune 
1698bf7f4e0aSPeter Brune #undef __FUNCT__
1699fa0ddf94SBarry Smith #define __FUNCT__ "SNESLineSearchSetWorkVecs"
17008d359177SBarry Smith /*@C
1701fa0ddf94SBarry Smith    SNESLineSearchSetWorkVecs - Gets work vectors for the line search.
1702f40b411bSPeter Brune 
1703f40b411bSPeter Brune    Input Parameter:
1704f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
1705f40b411bSPeter Brune -  nwork - the number of work vectors
1706f40b411bSPeter Brune 
1707f40b411bSPeter Brune    Level: developer
1708f40b411bSPeter Brune 
1709fa0ddf94SBarry Smith    Developers Note: This is PETSC_EXTERN because it may be used by user written plugin SNES implementations
1710cd7522eaSPeter Brune 
1711f1c6b773SPeter Brune .keywords: SNESLineSearch, work, vector
1712f40b411bSPeter Brune 
1713fa0ddf94SBarry Smith .seealso: SNESSetWorkVecs()
1714f40b411bSPeter Brune @*/
1715fa0ddf94SBarry Smith PetscErrorCode  SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork)
1716bf7f4e0aSPeter Brune {
1717bf7f4e0aSPeter Brune   PetscErrorCode ierr;
1718bf388a1fSBarry Smith 
1719bf7f4e0aSPeter Brune   PetscFunctionBegin;
1720bf7f4e0aSPeter Brune   if (linesearch->vec_sol) {
1721bf7f4e0aSPeter Brune     ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr);
17228d359177SBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
1723bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1724bf7f4e0aSPeter Brune }
1725bf7f4e0aSPeter Brune 
1726bf7f4e0aSPeter Brune #undef __FUNCT__
1727422a814eSBarry Smith #define __FUNCT__ "SNESLineSearchGetReason"
1728f40b411bSPeter Brune /*@
1729422a814eSBarry Smith    SNESLineSearchGetReason - Gets the success/failure status of the last line search application
1730f40b411bSPeter Brune 
1731f40b411bSPeter Brune    Input Parameters:
173278bcb3b5SPeter Brune .  linesearch - linesearch context
1733f40b411bSPeter Brune 
1734f40b411bSPeter Brune    Output Parameters:
1735422a814eSBarry Smith .  result - The success or failure status
1736f40b411bSPeter Brune 
1737cd7522eaSPeter Brune    Notes:
1738c98378a5SBarry Smith    This is typically called after SNESLineSearchApply() in order to determine if the line-search failed
1739cd7522eaSPeter Brune    (and set the SNES convergence accordingly).
1740cd7522eaSPeter Brune 
1741f40b411bSPeter Brune    Level: intermediate
1742f40b411bSPeter Brune 
1743422a814eSBarry Smith .seealso: SNESLineSearchSetReason(), SNESLineSearchReason
1744f40b411bSPeter Brune @*/
1745422a814eSBarry Smith PetscErrorCode  SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result)
1746bf7f4e0aSPeter Brune {
1747bf7f4e0aSPeter Brune   PetscFunctionBegin;
1748f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1749422a814eSBarry Smith   PetscValidPointer(result, 2);
1750422a814eSBarry Smith   *result = linesearch->result;
1751bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1752bf7f4e0aSPeter Brune }
1753bf7f4e0aSPeter Brune 
1754bf7f4e0aSPeter Brune #undef __FUNCT__
1755422a814eSBarry Smith #define __FUNCT__ "SNESLineSearchSetReason"
1756f40b411bSPeter Brune /*@
1757422a814eSBarry Smith    SNESLineSearchSetReason - Sets the success/failure status of the last line search application
1758f40b411bSPeter Brune 
1759f40b411bSPeter Brune    Input Parameters:
176078bcb3b5SPeter Brune +  linesearch - linesearch context
1761422a814eSBarry Smith -  result - The success or failure status
1762f40b411bSPeter Brune 
1763cd7522eaSPeter Brune    Notes:
1764cd7522eaSPeter Brune    This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set
1765cd7522eaSPeter Brune    the success or failure of the line search method.
1766cd7522eaSPeter Brune 
176778bcb3b5SPeter Brune    Level: developer
1768f40b411bSPeter Brune 
1769422a814eSBarry Smith .seealso: SNESLineSearchGetSResult()
1770f40b411bSPeter Brune @*/
1771422a814eSBarry Smith PetscErrorCode  SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result)
17726a388c36SPeter Brune {
17736a388c36SPeter Brune   PetscFunctionBegin;
17745fd66863SKarl Rupp   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1775422a814eSBarry Smith   linesearch->result = result;
17766a388c36SPeter Brune   PetscFunctionReturn(0);
17776a388c36SPeter Brune }
17786a388c36SPeter Brune 
17796a388c36SPeter Brune #undef __FUNCT__
1780f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetVIFunctions"
17819bd66eb0SPeter Brune /*@C
1782f1c6b773SPeter Brune    SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.
17839bd66eb0SPeter Brune 
17849bd66eb0SPeter Brune    Input Parameters:
17859bd66eb0SPeter Brune +  snes - nonlinear context obtained from SNESCreate()
17869bd66eb0SPeter Brune .  projectfunc - function for projecting the function to the bounds
17879bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
17889bd66eb0SPeter Brune 
17899bd66eb0SPeter Brune    Logically Collective on SNES
17909bd66eb0SPeter Brune 
17919bd66eb0SPeter Brune    Calling sequence of projectfunc:
17929bd66eb0SPeter Brune .vb
17939bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X)
17949bd66eb0SPeter Brune .ve
17959bd66eb0SPeter Brune 
17969bd66eb0SPeter Brune     Input parameters for projectfunc:
17979bd66eb0SPeter Brune +   snes - nonlinear context
17989bd66eb0SPeter Brune -   X - current solution
17999bd66eb0SPeter Brune 
1800cd7522eaSPeter Brune     Output parameters for projectfunc:
18019bd66eb0SPeter Brune .   X - Projected solution
18029bd66eb0SPeter Brune 
18039bd66eb0SPeter Brune    Calling sequence of normfunc:
18049bd66eb0SPeter Brune .vb
18059bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm)
18069bd66eb0SPeter Brune .ve
18079bd66eb0SPeter Brune 
1808cd7522eaSPeter Brune     Input parameters for normfunc:
18099bd66eb0SPeter Brune +   snes - nonlinear context
18109bd66eb0SPeter Brune .   X - current solution
18119bd66eb0SPeter Brune -   F - current residual
18129bd66eb0SPeter Brune 
1813cd7522eaSPeter Brune     Output parameters for normfunc:
18149bd66eb0SPeter Brune .   fnorm - VI-specific norm of the function
18159bd66eb0SPeter Brune 
1816cd7522eaSPeter Brune     Notes:
1817cd7522eaSPeter Brune     The VI solvers require projection of the solution to the feasible set.  projectfunc should implement this.
1818cd7522eaSPeter Brune 
1819cd7522eaSPeter Brune     The VI solvers require special evaluation of the function norm such that the norm is only calculated
1820cd7522eaSPeter Brune     on the inactive set.  This should be implemented by normfunc.
18219bd66eb0SPeter Brune 
18229bd66eb0SPeter Brune     Level: developer
18239bd66eb0SPeter Brune 
18249bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, set, line search
18259bd66eb0SPeter Brune 
1826f1c6b773SPeter Brune .seealso: SNESLineSearchGetVIFunctions(), SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck()
18279bd66eb0SPeter Brune @*/
1828f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc)
18299bd66eb0SPeter Brune {
18309bd66eb0SPeter Brune   PetscFunctionBegin;
1831f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
18329bd66eb0SPeter Brune   if (projectfunc) linesearch->ops->viproject = projectfunc;
18339bd66eb0SPeter Brune   if (normfunc) linesearch->ops->vinorm = normfunc;
18349bd66eb0SPeter Brune   PetscFunctionReturn(0);
18359bd66eb0SPeter Brune }
18369bd66eb0SPeter Brune 
18379bd66eb0SPeter Brune #undef __FUNCT__
1838f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetVIFunctions"
18399bd66eb0SPeter Brune /*@C
1840f1c6b773SPeter Brune    SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.
18419bd66eb0SPeter Brune 
18429bd66eb0SPeter Brune    Input Parameters:
1843907376e6SBarry Smith .  linesearch - the line search context, obtain with SNESGetLineSearch()
18449bd66eb0SPeter Brune 
18459bd66eb0SPeter Brune    Output Parameters:
18469bd66eb0SPeter Brune +  projectfunc - function for projecting the function to the bounds
18479bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
18489bd66eb0SPeter Brune 
18499bd66eb0SPeter Brune    Logically Collective on SNES
18509bd66eb0SPeter Brune 
18519bd66eb0SPeter Brune     Level: developer
18529bd66eb0SPeter Brune 
18539bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, get, line search
18549bd66eb0SPeter Brune 
1855f1c6b773SPeter Brune .seealso: SNESLineSearchSetVIFunctions(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
18569bd66eb0SPeter Brune @*/
1857f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc)
18589bd66eb0SPeter Brune {
18599bd66eb0SPeter Brune   PetscFunctionBegin;
18609bd66eb0SPeter Brune   if (projectfunc) *projectfunc = linesearch->ops->viproject;
18619bd66eb0SPeter Brune   if (normfunc) *normfunc = linesearch->ops->vinorm;
18629bd66eb0SPeter Brune   PetscFunctionReturn(0);
18639bd66eb0SPeter Brune }
18649bd66eb0SPeter Brune 
18659bd66eb0SPeter Brune #undef __FUNCT__
1866f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchRegister"
1867bf7f4e0aSPeter Brune /*@C
18681c84c290SBarry Smith   SNESLineSearchRegister - See SNESLineSearchRegister()
1869bf7f4e0aSPeter Brune 
1870bf7f4e0aSPeter Brune   Level: advanced
1871bf7f4e0aSPeter Brune @*/
1872bdf89e91SBarry Smith PetscErrorCode  SNESLineSearchRegister(const char sname[],PetscErrorCode (*function)(SNESLineSearch))
1873bf7f4e0aSPeter Brune {
1874bf7f4e0aSPeter Brune   PetscErrorCode ierr;
1875bf7f4e0aSPeter Brune 
1876bf7f4e0aSPeter Brune   PetscFunctionBegin;
1877a240a19fSJed Brown   ierr = PetscFunctionListAdd(&SNESLineSearchList,sname,function);CHKERRQ(ierr);
1878bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1879bf7f4e0aSPeter Brune }
1880