xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision cd7522eaf66b5dfe973598b9eb081fb87c780851)
19e764e56SPeter Brune #include <private/linesearchimpl.h> /*I "petscsnes.h" I*/
2bf7f4e0aSPeter Brune 
3f1c6b773SPeter Brune PetscBool  SNESLineSearchRegisterAllCalled = PETSC_FALSE;
4f1c6b773SPeter Brune PetscFList SNESLineSearchList              = PETSC_NULL;
5bf7f4e0aSPeter Brune 
6f1c6b773SPeter Brune PetscClassId   SNESLINESEARCH_CLASSID;
7f1c6b773SPeter Brune PetscLogEvent  SNESLineSearch_Apply;
8bf7f4e0aSPeter Brune 
9bf7f4e0aSPeter Brune #undef __FUNCT__
10f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate"
11f40b411bSPeter Brune /*@
12*cd7522eaSPeter Brune    SNESLineSearchCreate - Creates the line search context.
13f40b411bSPeter Brune 
14*cd7522eaSPeter Brune    Logically Collective on Comm
15f40b411bSPeter Brune 
16f40b411bSPeter Brune    Input Parameters:
17*cd7522eaSPeter Brune .  comm - MPI communicator for the line search (typically from the associated SNES context).
18f40b411bSPeter Brune 
19f40b411bSPeter Brune    Output Parameters:
20f40b411bSPeter Brune .  outlinesearch - the line search instance.
21f40b411bSPeter Brune 
22f40b411bSPeter Brune    Level: Beginner
23f40b411bSPeter Brune 
24*cd7522eaSPeter Brune    .keywords: LineSearch, create, context
25f40b411bSPeter Brune 
26f40b411bSPeter Brune    .seealso: LineSearchDestroy()
27f40b411bSPeter Brune @*/
28f40b411bSPeter Brune 
29f1c6b773SPeter Brune PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch) {
30bf7f4e0aSPeter Brune   PetscErrorCode      ierr;
31f1c6b773SPeter Brune   SNESLineSearch     linesearch;
32bf7f4e0aSPeter Brune   PetscFunctionBegin;
33ea5d4fccSPeter Brune   PetscValidPointer(outlinesearch,2);
34ea5d4fccSPeter Brune   *outlinesearch = PETSC_NULL;
35f1c6b773SPeter Brune   ierr = PetscHeaderCreate(linesearch,_p_LineSearch,struct _LineSearchOps,SNESLINESEARCH_CLASSID, 0,
36f1c6b773SPeter Brune                            "SNESLineSearch","Line-search method","SNESLineSearch",comm,SNESLineSearchDestroy,SNESLineSearchView);CHKERRQ(ierr);
37bf7f4e0aSPeter Brune 
38bf7f4e0aSPeter Brune   linesearch->ops->precheckstep = PETSC_NULL;
39bf7f4e0aSPeter Brune   linesearch->ops->postcheckstep = PETSC_NULL;
40bf7f4e0aSPeter Brune 
419bd66eb0SPeter Brune   linesearch->vec_sol_new   = PETSC_NULL;
429bd66eb0SPeter Brune   linesearch->vec_func_new  = PETSC_NULL;
439bd66eb0SPeter Brune   linesearch->vec_sol       = PETSC_NULL;
449bd66eb0SPeter Brune   linesearch->vec_func      = PETSC_NULL;
459bd66eb0SPeter Brune   linesearch->vec_update    = PETSC_NULL;
469bd66eb0SPeter Brune 
47bf7f4e0aSPeter Brune   linesearch->lambda        = 1.0;
48bf7f4e0aSPeter Brune   linesearch->fnorm         = 1.0;
49bf7f4e0aSPeter Brune   linesearch->ynorm         = 1.0;
50bf7f4e0aSPeter Brune   linesearch->xnorm         = 1.0;
51bf7f4e0aSPeter Brune   linesearch->success       = PETSC_TRUE;
52bf7f4e0aSPeter Brune   linesearch->norms         = PETSC_TRUE;
53bf7f4e0aSPeter Brune   linesearch->keeplambda    = PETSC_FALSE;
54bf7f4e0aSPeter Brune   linesearch->damping       = 1.0;
55bf7f4e0aSPeter Brune   linesearch->maxstep       = 1e8;
56bf7f4e0aSPeter Brune   linesearch->steptol       = 1e-12;
57516fe3c3SPeter Brune   linesearch->rtol          = 1e-8;
58516fe3c3SPeter Brune   linesearch->atol          = 1e-15;
59516fe3c3SPeter Brune   linesearch->ltol          = 1e-8;
60bf7f4e0aSPeter Brune   linesearch->precheckctx   = PETSC_NULL;
61bf7f4e0aSPeter Brune   linesearch->postcheckctx  = PETSC_NULL;
62516fe3c3SPeter Brune   linesearch->max_its       = 3;
63bf7f4e0aSPeter Brune   linesearch->setupcalled   = PETSC_FALSE;
64bf7f4e0aSPeter Brune   *outlinesearch            = linesearch;
65bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
66bf7f4e0aSPeter Brune }
67bf7f4e0aSPeter Brune 
68bf7f4e0aSPeter Brune #undef __FUNCT__
69f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetUp"
70f40b411bSPeter Brune /*@
71f1c6b773SPeter Brune    SNESLineSearchSetUp - Prepares the line search for being applied.
72f40b411bSPeter Brune 
73*cd7522eaSPeter Brune    Collective on SNESLineSearch
74f40b411bSPeter Brune 
75f40b411bSPeter Brune    Input Parameters:
76f40b411bSPeter Brune .  linesearch - The LineSearch instance.
77f40b411bSPeter Brune 
78*cd7522eaSPeter Brune    Notes:
79*cd7522eaSPeter Brune    For most cases, this needn't be called outside of SNESLineSearchApply().
80*cd7522eaSPeter Brune    The only current case where this is called outside of this is for the VI
81*cd7522eaSPeter Brune    solvers, which modifies the solution and work vectors before the first call
82*cd7522eaSPeter Brune    of SNESLineSearchApply, requiring the SNESLineSearch work vectors to be
83*cd7522eaSPeter Brune    allocated upfront.
84*cd7522eaSPeter Brune 
85*cd7522eaSPeter Brune 
86f40b411bSPeter Brune    Level: Intermediate
87f40b411bSPeter Brune 
88f1c6b773SPeter Brune    .keywords: SNESLineSearch, SetUp
89f40b411bSPeter Brune 
90f1c6b773SPeter Brune    .seealso: SNESLineSearchReset()
91f40b411bSPeter Brune @*/
92f40b411bSPeter Brune 
93f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch) {
94bf7f4e0aSPeter Brune   PetscErrorCode ierr;
95bf7f4e0aSPeter Brune   PetscFunctionBegin;
96bf7f4e0aSPeter Brune   if (!((PetscObject)linesearch)->type_name) {
97f1c6b773SPeter Brune     ierr = SNESLineSearchSetType(linesearch,SNES_LINESEARCH_BASIC);CHKERRQ(ierr);
98bf7f4e0aSPeter Brune   }
99bf7f4e0aSPeter Brune   if (!linesearch->setupcalled) {
1009bd66eb0SPeter Brune     if (!linesearch->vec_sol_new) {
101bf7f4e0aSPeter Brune       ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new);CHKERRQ(ierr);
1029bd66eb0SPeter Brune     }
1039bd66eb0SPeter Brune     if (!linesearch->vec_func_new) {
1049bd66eb0SPeter Brune       ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new);CHKERRQ(ierr);
1059bd66eb0SPeter Brune     }
106bf7f4e0aSPeter Brune     if (linesearch->ops->setup) {
107bf7f4e0aSPeter Brune       ierr = (*linesearch->ops->setup)(linesearch);CHKERRQ(ierr);
108bf7f4e0aSPeter Brune     }
109bf7f4e0aSPeter Brune     linesearch->lambda = linesearch->damping;
110bf7f4e0aSPeter Brune     linesearch->setupcalled = PETSC_TRUE;
111bf7f4e0aSPeter Brune   }
112bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
113bf7f4e0aSPeter Brune }
114bf7f4e0aSPeter Brune 
115bf7f4e0aSPeter Brune #undef __FUNCT__
116f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchReset"
117f40b411bSPeter Brune 
118f40b411bSPeter Brune /*@
119*cd7522eaSPeter Brune    SNESLineSearchReset - Tears down the structures required for application of the SNESLineSearch
120f40b411bSPeter Brune 
121f1c6b773SPeter Brune    Collective on SNESLineSearch
122f40b411bSPeter Brune 
123f40b411bSPeter Brune    Input Parameters:
124f40b411bSPeter Brune .  linesearch - The LineSearch instance.
125f40b411bSPeter Brune 
126f40b411bSPeter Brune    Level: Intermediate
127f40b411bSPeter Brune 
128*cd7522eaSPeter Brune    .keywords: SNESLineSearch, Reset
129f40b411bSPeter Brune 
130f1c6b773SPeter Brune    .seealso: SNESLineSearchSetUp()
131f40b411bSPeter Brune @*/
132f40b411bSPeter Brune 
133f1c6b773SPeter Brune PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch) {
134bf7f4e0aSPeter Brune   PetscErrorCode ierr;
135bf7f4e0aSPeter Brune   PetscFunctionBegin;
136bf7f4e0aSPeter Brune   if (linesearch->ops->reset) {
137bf7f4e0aSPeter Brune     (*linesearch->ops->reset)(linesearch);
138bf7f4e0aSPeter Brune   }
139bf7f4e0aSPeter Brune   ierr = VecDestroy(&linesearch->vec_sol_new);CHKERRQ(ierr);
140bf7f4e0aSPeter Brune   ierr = VecDestroy(&linesearch->vec_func_new);CHKERRQ(ierr);
141bf7f4e0aSPeter Brune 
142bf7f4e0aSPeter Brune   ierr = VecDestroyVecs(linesearch->nwork, &linesearch->work);CHKERRQ(ierr);
143bf7f4e0aSPeter Brune   linesearch->nwork = 0;
144bf7f4e0aSPeter Brune   linesearch->setupcalled = PETSC_FALSE;
145bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
146bf7f4e0aSPeter Brune }
147bf7f4e0aSPeter Brune 
14886d74e61SPeter Brune 
14986d74e61SPeter Brune #undef __FUNCT__
150f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetPreCheck"
15186d74e61SPeter Brune /*@C
152f1c6b773SPeter Brune    SNESLineSearchSetPreCheck - Sets a pre-check function for the line search routine.
15386d74e61SPeter Brune 
154f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
15586d74e61SPeter Brune 
15686d74e61SPeter Brune    Input Parameters:
157f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
15886d74e61SPeter Brune .  func       - [optional] function evaluation routine
15986d74e61SPeter Brune -  ctx        - [optional] user-defined context for private data for the
16086d74e61SPeter Brune                 function evaluation routine (may be PETSC_NULL)
16186d74e61SPeter Brune 
16286d74e61SPeter Brune    Calling sequence of func:
163f1c6b773SPeter Brune $    func (SNESLineSearch snes,Vec x,Vec y, PetscBool *changed);
16486d74e61SPeter Brune 
16586d74e61SPeter Brune +  x - solution vector
16686d74e61SPeter Brune .  y - search direction vector
167*cd7522eaSPeter Brune -  changed - flag to indicate the precheck changed x or y.
16886d74e61SPeter Brune 
16986d74e61SPeter Brune    Level: intermediate
17086d74e61SPeter Brune 
17186d74e61SPeter Brune .keywords: set, linesearch, pre-check
17286d74e61SPeter Brune 
173f1c6b773SPeter Brune .seealso: SNESLineSearchSetPostCheck()
17486d74e61SPeter Brune @*/
175f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetPreCheck(SNESLineSearch linesearch, SNESLineSearchPreCheckFunc func,void *ctx)
17686d74e61SPeter Brune {
1779bd66eb0SPeter Brune   PetscFunctionBegin;
178f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
17986d74e61SPeter Brune   if (func) linesearch->ops->precheckstep = func;
18086d74e61SPeter Brune   if (ctx) linesearch->precheckctx = ctx;
18186d74e61SPeter Brune   PetscFunctionReturn(0);
18286d74e61SPeter Brune }
18386d74e61SPeter Brune 
18486d74e61SPeter Brune 
18586d74e61SPeter Brune #undef __FUNCT__
186f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetPreCheck"
18786d74e61SPeter Brune /*@C
188*cd7522eaSPeter Brune    SNESLineSearchSetPreCheck - Gets the pre-check function for the line search routine.
18986d74e61SPeter Brune 
19086d74e61SPeter Brune    Input Parameters:
191f1c6b773SPeter Brune .  linesearch - the SNESLineSearch context
19286d74e61SPeter Brune 
19386d74e61SPeter Brune    Output Parameters:
19486d74e61SPeter Brune +  func       - [optional] function evaluation routine
19586d74e61SPeter Brune -  ctx        - [optional] user-defined context for private data for the
19686d74e61SPeter Brune                 function evaluation routine (may be PETSC_NULL)
19786d74e61SPeter Brune 
19886d74e61SPeter Brune    Level: intermediate
19986d74e61SPeter Brune 
20086d74e61SPeter Brune .keywords: get, linesearch, pre-check
20186d74e61SPeter Brune 
202f1c6b773SPeter Brune .seealso: SNESLineSearchGetPostCheck(), SNESLineSearchSetPreCheck()
20386d74e61SPeter Brune @*/
204f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetPreCheck(SNESLineSearch linesearch, SNESLineSearchPreCheckFunc *func,void **ctx)
20586d74e61SPeter Brune {
2069bd66eb0SPeter Brune   PetscFunctionBegin;
207f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
20886d74e61SPeter Brune   if (func) *func = linesearch->ops->precheckstep;
20986d74e61SPeter Brune   if (ctx) *ctx = linesearch->precheckctx;
21086d74e61SPeter Brune   PetscFunctionReturn(0);
21186d74e61SPeter Brune }
21286d74e61SPeter Brune 
21386d74e61SPeter Brune 
21486d74e61SPeter Brune #undef __FUNCT__
215f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetPostCheck"
21686d74e61SPeter Brune /*@C
217f1c6b773SPeter Brune    SNESLineSearchSetPostCheck - Sets a post-check function for the line search routine.
21886d74e61SPeter Brune 
219f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
22086d74e61SPeter Brune 
22186d74e61SPeter Brune    Input Parameters:
222f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
22386d74e61SPeter Brune .  func       - [optional] function evaluation routine
22486d74e61SPeter Brune -  ctx        - [optional] user-defined context for private data for the
22586d74e61SPeter Brune                 function evaluation routine (may be PETSC_NULL)
22686d74e61SPeter Brune 
22786d74e61SPeter Brune    Calling sequence of func:
228f1c6b773SPeter Brune $    func (SNESLineSearch linesearch,Vec x, Vec w, Vec y, PetscBool *changed_w, *changed_y);
22986d74e61SPeter Brune 
23086d74e61SPeter Brune +  x - old solution vector
23186d74e61SPeter Brune .  y - search direction vector
23286d74e61SPeter Brune .  w - new solution vector
23386d74e61SPeter Brune .  changed_y - indicates that the line search changed y.
23486d74e61SPeter Brune .  changed_w - indicates that the line search changed w.
23586d74e61SPeter Brune 
23686d74e61SPeter Brune    Level: intermediate
23786d74e61SPeter Brune 
23886d74e61SPeter Brune .keywords: set, linesearch, post-check
23986d74e61SPeter Brune 
240f1c6b773SPeter Brune .seealso: SNESLineSearchSetPreCheck()
24186d74e61SPeter Brune @*/
242f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetPostCheck(SNESLineSearch linesearch, SNESLineSearchPostCheckFunc func,void *ctx)
24386d74e61SPeter Brune {
24486d74e61SPeter Brune   PetscFunctionBegin;
245f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
24686d74e61SPeter Brune   if (func) linesearch->ops->postcheckstep = func;
24786d74e61SPeter Brune   if (ctx) linesearch->postcheckctx = ctx;
24886d74e61SPeter Brune   PetscFunctionReturn(0);
24986d74e61SPeter Brune }
25086d74e61SPeter Brune 
25186d74e61SPeter Brune 
25286d74e61SPeter Brune #undef __FUNCT__
253f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetPostCheck"
25486d74e61SPeter Brune /*@C
255f1c6b773SPeter Brune    SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine.
25686d74e61SPeter Brune 
25786d74e61SPeter Brune    Input Parameters:
258f1c6b773SPeter Brune .  linesearch - the SNESLineSearch context
25986d74e61SPeter Brune 
26086d74e61SPeter Brune    Output Parameters:
26186d74e61SPeter Brune +  func       - [optional] function evaluation routine
26286d74e61SPeter Brune -  ctx        - [optional] user-defined context for private data for the
26386d74e61SPeter Brune                 function evaluation routine (may be PETSC_NULL)
26486d74e61SPeter Brune 
26586d74e61SPeter Brune    Level: intermediate
26686d74e61SPeter Brune 
26786d74e61SPeter Brune .keywords: get, linesearch, post-check
26886d74e61SPeter Brune 
269f1c6b773SPeter Brune .seealso: SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck()
27086d74e61SPeter Brune @*/
271f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetPostCheck(SNESLineSearch linesearch, SNESLineSearchPostCheckFunc *func,void **ctx)
27286d74e61SPeter Brune {
2739bd66eb0SPeter Brune   PetscFunctionBegin;
274f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
27586d74e61SPeter Brune   if (func) *func = linesearch->ops->postcheckstep;
27686d74e61SPeter Brune   if (ctx) *ctx = linesearch->postcheckctx;
27786d74e61SPeter Brune   PetscFunctionReturn(0);
27886d74e61SPeter Brune }
27986d74e61SPeter Brune 
28086d74e61SPeter Brune 
281bf7f4e0aSPeter Brune #undef __FUNCT__
282f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPreCheck"
283f40b411bSPeter Brune /*@
284f1c6b773SPeter Brune    SNESLineSearchPreCheck - Prepares the line search for being applied.
285f40b411bSPeter Brune 
286*cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
287f40b411bSPeter Brune 
288f40b411bSPeter Brune    Input Parameters:
289f40b411bSPeter Brune .  linesearch - The linesearch instance.
290f40b411bSPeter Brune 
291f40b411bSPeter Brune    Output Parameters:
292f40b411bSPeter Brune .  changed - Indicator if the pre-check has changed anything.
293f40b411bSPeter Brune 
294f40b411bSPeter Brune    Level: Beginner
295f40b411bSPeter Brune 
296f1c6b773SPeter Brune    .keywords: SNESLineSearch, Create
297f40b411bSPeter Brune 
298f1c6b773SPeter Brune    .seealso: SNESLineSearchPostCheck()
299f40b411bSPeter Brune @*/
300f1c6b773SPeter Brune PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch, PetscBool * changed)
301bf7f4e0aSPeter Brune {
302bf7f4e0aSPeter Brune   PetscErrorCode ierr;
303bf7f4e0aSPeter Brune   PetscFunctionBegin;
304bf7f4e0aSPeter Brune   *changed = PETSC_FALSE;
305bf7f4e0aSPeter Brune   if (linesearch->ops->precheckstep) {
306c87759e9SPeter Brune     ierr = (*linesearch->ops->precheckstep)(linesearch, linesearch->vec_sol, linesearch->vec_update, changed, linesearch->precheckctx);CHKERRQ(ierr);
307bf7f4e0aSPeter Brune   }
308bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
309bf7f4e0aSPeter Brune }
310bf7f4e0aSPeter Brune 
311bf7f4e0aSPeter Brune #undef __FUNCT__
312f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPostCheck"
313f40b411bSPeter Brune /*@
314f1c6b773SPeter Brune    SNESLineSearchPostCheck - Prepares the line search for being applied.
315f40b411bSPeter Brune 
316*cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
317f40b411bSPeter Brune 
318f40b411bSPeter Brune    Input Parameters:
319f40b411bSPeter Brune .  linesearch - The linesearch instance.
320f40b411bSPeter Brune 
321f40b411bSPeter Brune    Output Parameters:
32268fb48e8SJed Brown +  changed_Y - Indicator if the direction has been changed.
32368fb48e8SJed Brown -  changed_W - Indicator if the solution has been changed.
324f40b411bSPeter Brune 
325f40b411bSPeter Brune    Level: Intermediate
326f40b411bSPeter Brune 
327f1c6b773SPeter Brune    .keywords: SNESLineSearch, Create
328f40b411bSPeter Brune 
329f1c6b773SPeter Brune    .seealso: SNESLineSearchPreCheck()
330f40b411bSPeter Brune @*/
331f1c6b773SPeter Brune PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch, PetscBool * changed_Y, PetscBool * changed_W)
332bf7f4e0aSPeter Brune {
333bf7f4e0aSPeter Brune   PetscErrorCode ierr;
334bf7f4e0aSPeter Brune   PetscFunctionBegin;
335bf7f4e0aSPeter Brune   *changed_Y = PETSC_FALSE;
336bf7f4e0aSPeter Brune   *changed_W = PETSC_FALSE;
337bf7f4e0aSPeter Brune   if (linesearch->ops->postcheckstep) {
338c87759e9SPeter Brune     ierr = (*linesearch->ops->postcheckstep)(linesearch, linesearch->vec_sol, linesearch->vec_update, linesearch->vec_sol_new, changed_Y, changed_W, linesearch->postcheckctx);CHKERRQ(ierr);
33986d74e61SPeter Brune   }
34086d74e61SPeter Brune   PetscFunctionReturn(0);
34186d74e61SPeter Brune }
34286d74e61SPeter Brune 
34386d74e61SPeter Brune 
34486d74e61SPeter Brune #undef __FUNCT__
345f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPreCheckPicard"
34686d74e61SPeter Brune /*@C
34786d74e61SPeter Brune    SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration
34886d74e61SPeter Brune 
349*cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
35086d74e61SPeter Brune 
35186d74e61SPeter Brune    Input Arguments:
35286d74e61SPeter Brune +  linesearch - linesearch context
35386d74e61SPeter Brune .  X - base state for this step
35486d74e61SPeter Brune .  Y - initial correction
35586d74e61SPeter Brune 
35686d74e61SPeter Brune    Output Arguments:
35786d74e61SPeter Brune +  Y - correction, possibly modified
35886d74e61SPeter Brune -  changed - flag indicating that Y was modified
35986d74e61SPeter Brune 
36086d74e61SPeter Brune    Options Database Key:
361*cd7522eaSPeter Brune +  -snes_linesearch_precheck_picard - activate this routine
362*cd7522eaSPeter Brune -  -snes_linesearch_precheck_picard_angle - angle
36386d74e61SPeter Brune 
36486d74e61SPeter Brune    Level: advanced
36586d74e61SPeter Brune 
36686d74e61SPeter Brune    Notes:
36786d74e61SPeter Brune    This function should be passed to SNESLineSearchSetPreCheck()
36886d74e61SPeter Brune 
36986d74e61SPeter Brune    The justification for this method involves the linear convergence of a Picard iteration
37086d74e61SPeter Brune    so the Picard linearization should be provided in place of the "Jacobian". This correction
37186d74e61SPeter Brune    is generally not useful when using a Newton linearization.
37286d74e61SPeter Brune 
37386d74e61SPeter Brune    Reference:
37486d74e61SPeter Brune    Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology.
37586d74e61SPeter Brune 
37686d74e61SPeter Brune .seealso: SNESLineSearchSetPreCheck()
37786d74e61SPeter Brune @*/
378f1c6b773SPeter Brune PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed,void *ctx)
37986d74e61SPeter Brune {
38086d74e61SPeter Brune   PetscErrorCode ierr;
38186d74e61SPeter Brune   PetscReal      angle = *(PetscReal*)linesearch->precheckctx;
38286d74e61SPeter Brune   Vec            Ylast;
38386d74e61SPeter Brune   PetscScalar    dot;
384c87759e9SPeter Brune 
38586d74e61SPeter Brune   PetscInt       iter;
38686d74e61SPeter Brune   PetscReal      ynorm,ylastnorm,theta,angle_radians;
38786d74e61SPeter Brune   SNES           snes;
38886d74e61SPeter Brune 
38986d74e61SPeter Brune   PetscFunctionBegin;
390f1c6b773SPeter Brune   ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
39186d74e61SPeter Brune   ierr = PetscObjectQuery((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject*)&Ylast);CHKERRQ(ierr);
39286d74e61SPeter Brune   if (!Ylast) {
39386d74e61SPeter Brune     ierr = VecDuplicate(Y,&Ylast);CHKERRQ(ierr);
39486d74e61SPeter Brune     ierr = PetscObjectCompose((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject)Ylast);CHKERRQ(ierr);
39586d74e61SPeter Brune     ierr = PetscObjectDereference((PetscObject)Ylast);CHKERRQ(ierr);
39686d74e61SPeter Brune   }
39786d74e61SPeter Brune   ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr);
39886d74e61SPeter Brune   if (iter < 2) {
39986d74e61SPeter Brune     ierr = VecCopy(Y,Ylast);CHKERRQ(ierr);
40086d74e61SPeter Brune     *changed = PETSC_FALSE;
40186d74e61SPeter Brune     PetscFunctionReturn(0);
40286d74e61SPeter Brune   }
40386d74e61SPeter Brune 
40486d74e61SPeter Brune   ierr = VecDot(Y,Ylast,&dot);CHKERRQ(ierr);
40586d74e61SPeter Brune   ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr);
40686d74e61SPeter Brune   ierr = VecNorm(Ylast,NORM_2,&ylastnorm);CHKERRQ(ierr);
40786d74e61SPeter Brune   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
40886d74e61SPeter Brune   theta = acos((double)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0));
40986d74e61SPeter Brune   angle_radians = angle * PETSC_PI / 180.;
41086d74e61SPeter Brune   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
41186d74e61SPeter Brune     /* Modify the step Y */
41286d74e61SPeter Brune     PetscReal alpha,ydiffnorm;
41386d74e61SPeter Brune     ierr = VecAXPY(Ylast,-1.0,Y);CHKERRQ(ierr);
41486d74e61SPeter Brune     ierr = VecNorm(Ylast,NORM_2,&ydiffnorm);CHKERRQ(ierr);
41586d74e61SPeter Brune     alpha = ylastnorm / ydiffnorm;
41686d74e61SPeter Brune     ierr = VecCopy(Y,Ylast);CHKERRQ(ierr);
41786d74e61SPeter Brune     ierr = VecScale(Y,alpha);CHKERRQ(ierr);
41886d74e61SPeter Brune     ierr = PetscInfo3(snes,"Angle %G degrees less than threshold %G, corrected step by alpha=%G\n",theta*180./PETSC_PI,angle,alpha);CHKERRQ(ierr);
41986d74e61SPeter Brune   } else {
42086d74e61SPeter Brune     ierr = PetscInfo2(snes,"Angle %G degrees exceeds threshold %G, no correction applied\n",theta*180./PETSC_PI,angle);CHKERRQ(ierr);
42186d74e61SPeter Brune     ierr = VecCopy(Y,Ylast);CHKERRQ(ierr);
42286d74e61SPeter Brune     *changed = PETSC_FALSE;
423bf7f4e0aSPeter Brune   }
424bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
425bf7f4e0aSPeter Brune }
426bf7f4e0aSPeter Brune 
427bf7f4e0aSPeter Brune #undef __FUNCT__
428f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchApply"
429f40b411bSPeter Brune /*@
430*cd7522eaSPeter Brune    SNESLineSearchApply - Computes the line-search update.
431f40b411bSPeter Brune 
432f1c6b773SPeter Brune    Collective on SNESLineSearch
433f40b411bSPeter Brune 
434f40b411bSPeter Brune    Input Parameters:
435f40b411bSPeter Brune +  linesearch - The linesearch instance.
436f40b411bSPeter Brune .  X - The current solution.
437f40b411bSPeter Brune .  F - The current function.
438f40b411bSPeter Brune .  fnorm - The current norm.
439f40b411bSPeter Brune .  Y - The search direction.
440f40b411bSPeter Brune 
441f40b411bSPeter Brune    Output Parameters:
442f40b411bSPeter Brune +  X - The new solution.
443f40b411bSPeter Brune .  F - The new function.
444f40b411bSPeter Brune -  fnorm - The new function norm.
445f40b411bSPeter Brune 
446*cd7522eaSPeter Brune    Options Database Keys:
447*cd7522eaSPeter Brune + -snes_linesearch_type - The Line search method
448*cd7522eaSPeter Brune . -snes_linesearch_monitor - Print progress of line searches
449*cd7522eaSPeter Brune . -snes_linesearch_damping - The linesearch damping parameter.
450*cd7522eaSPeter Brune . -snes_linesearch_norms   - Turn on/off the linesearch norms
451*cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess.
452*cd7522eaSPeter Brune - -snes_linesearch_max_it - The number of iterations for iterative line searches.
453*cd7522eaSPeter Brune 
454*cd7522eaSPeter Brune    Notes:
455*cd7522eaSPeter Brune    This is typically called from within a SNESSolve() implementation in order to
456*cd7522eaSPeter Brune    help with convergence of the nonlinear method.  Various SNES types use line searches
457*cd7522eaSPeter Brune    in different ways, but the overarching theme is that a line search is used to determine
458*cd7522eaSPeter Brune    an optimal damping parameter of a step at each iteration of the method.  Each
459*cd7522eaSPeter Brune    application of the line search may invoke SNESComputeFunction several times, and
460*cd7522eaSPeter Brune    therefore may be fairly expensive.
461*cd7522eaSPeter Brune 
462f40b411bSPeter Brune    Level: Intermediate
463f40b411bSPeter Brune 
464f1c6b773SPeter Brune    .keywords: SNESLineSearch, Create
465f40b411bSPeter Brune 
466*cd7522eaSPeter Brune    .seealso: SNESLineSearchCreate(), SNESLineSearchPreCheck(), SNESLineSearchPostCheck(), SNESSolve(), SNESComputeFunction()
467f40b411bSPeter Brune @*/
468f1c6b773SPeter Brune PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y) {
469bf7f4e0aSPeter Brune   PetscErrorCode ierr;
470bf7f4e0aSPeter Brune   PetscFunctionBegin;
471bf7f4e0aSPeter Brune 
472bf7f4e0aSPeter Brune   /* check the pointers */
473f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
474bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
475bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
476bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
477bf7f4e0aSPeter Brune 
478bf7f4e0aSPeter Brune   linesearch->success = PETSC_TRUE;
479bf7f4e0aSPeter Brune 
480bf7f4e0aSPeter Brune   linesearch->vec_sol = X;
481bf7f4e0aSPeter Brune   linesearch->vec_update = Y;
482bf7f4e0aSPeter Brune   linesearch->vec_func = F;
483bf7f4e0aSPeter Brune 
484f1c6b773SPeter Brune   ierr = SNESLineSearchSetUp(linesearch);CHKERRQ(ierr);
485bf7f4e0aSPeter Brune 
486bf7f4e0aSPeter Brune   if (!linesearch->keeplambda)
487bf7f4e0aSPeter Brune     linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */
488bf7f4e0aSPeter Brune 
489bf7f4e0aSPeter Brune   if (fnorm) {
490bf7f4e0aSPeter Brune     linesearch->fnorm = *fnorm;
491bf7f4e0aSPeter Brune   } else {
492bf7f4e0aSPeter Brune     ierr = VecNorm(F, NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
493bf7f4e0aSPeter Brune   }
494bf7f4e0aSPeter Brune 
495f1c6b773SPeter Brune   ierr = PetscLogEventBegin(SNESLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
496bf7f4e0aSPeter Brune 
497bf7f4e0aSPeter Brune   ierr = (*linesearch->ops->apply)(linesearch);CHKERRQ(ierr);
498bf7f4e0aSPeter Brune 
499f1c6b773SPeter Brune   ierr = PetscLogEventEnd(SNESLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
500bf7f4e0aSPeter Brune 
501bf7f4e0aSPeter Brune   if (fnorm)
502bf7f4e0aSPeter Brune     *fnorm = linesearch->fnorm;
503bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
504bf7f4e0aSPeter Brune }
505bf7f4e0aSPeter Brune 
506bf7f4e0aSPeter Brune #undef __FUNCT__
507f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchDestroy"
508f40b411bSPeter Brune /*@
509f1c6b773SPeter Brune    SNESLineSearchDestroy - Destroys the line search instance.
510f40b411bSPeter Brune 
511f1c6b773SPeter Brune    Collective on SNESLineSearch
512f40b411bSPeter Brune 
513f40b411bSPeter Brune    Input Parameters:
514f40b411bSPeter Brune .  linesearch - The linesearch instance.
515f40b411bSPeter Brune 
516f40b411bSPeter Brune    Level: Intermediate
517f40b411bSPeter Brune 
518f1c6b773SPeter Brune    .keywords: SNESLineSearch, Create
519f40b411bSPeter Brune 
520*cd7522eaSPeter Brune    .seealso: SNESLineSearchCreate(), SNESLineSearchReset(), SNESDestroy()
521f40b411bSPeter Brune @*/
522f1c6b773SPeter Brune PetscErrorCode SNESLineSearchDestroy(SNESLineSearch * linesearch) {
523bf7f4e0aSPeter Brune   PetscErrorCode ierr;
524bf7f4e0aSPeter Brune   PetscFunctionBegin;
525bf7f4e0aSPeter Brune   if (!*linesearch) PetscFunctionReturn(0);
526f1c6b773SPeter Brune   PetscValidHeaderSpecific((*linesearch),SNESLINESEARCH_CLASSID,1);
527bf7f4e0aSPeter Brune   if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = 0; PetscFunctionReturn(0);}
528bf7f4e0aSPeter Brune   ierr = PetscObjectDepublish((*linesearch));CHKERRQ(ierr);
529f1c6b773SPeter Brune   ierr = SNESLineSearchReset(*linesearch);
530bf7f4e0aSPeter Brune   if ((*linesearch)->ops->destroy) {
531bf7f4e0aSPeter Brune     (*linesearch)->ops->destroy(*linesearch);
532bf7f4e0aSPeter Brune   }
533bf7f4e0aSPeter Brune   ierr = PetscViewerDestroy(&(*linesearch)->monitor);CHKERRQ(ierr);
534e7058c64SPeter Brune   ierr = PetscHeaderDestroy(linesearch);CHKERRQ(ierr);
535bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
536bf7f4e0aSPeter Brune }
537bf7f4e0aSPeter Brune 
538bf7f4e0aSPeter Brune #undef __FUNCT__
539f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetMonitor"
540f40b411bSPeter Brune /*@
541*cd7522eaSPeter Brune    SNESLineSearchSetMonitor - Turns on/off printing useful information and debugging output about the line search.
542bf7f4e0aSPeter Brune 
543bf7f4e0aSPeter Brune    Input Parameters:
544bf7f4e0aSPeter Brune +  snes - nonlinear context obtained from SNESCreate()
545bf7f4e0aSPeter Brune -  flg - PETSC_TRUE to monitor the line search
546bf7f4e0aSPeter Brune 
547bf7f4e0aSPeter Brune    Logically Collective on SNES
548bf7f4e0aSPeter Brune 
549bf7f4e0aSPeter Brune    Options Database:
550*cd7522eaSPeter Brune .   -snes_linesearch_monitor - enables the monitor.
551bf7f4e0aSPeter Brune 
552bf7f4e0aSPeter Brune    Level: intermediate
553bf7f4e0aSPeter Brune 
554bf7f4e0aSPeter Brune 
555*cd7522eaSPeter Brune .seealso: SNESLineSearchGetMonitor(), PetscViewer
556bf7f4e0aSPeter Brune @*/
557f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetMonitor(SNESLineSearch linesearch, PetscBool flg)
558bf7f4e0aSPeter Brune {
559bf7f4e0aSPeter Brune 
560bf7f4e0aSPeter Brune   PetscErrorCode ierr;
561bf7f4e0aSPeter Brune   PetscFunctionBegin;
562bf7f4e0aSPeter Brune   if (flg && !linesearch->monitor) {
563bf7f4e0aSPeter Brune     ierr = PetscViewerASCIIOpen(((PetscObject)linesearch)->comm,"stdout",&linesearch->monitor);CHKERRQ(ierr);
564bf7f4e0aSPeter Brune   } else if (!flg && linesearch->monitor) {
565bf7f4e0aSPeter Brune     ierr = PetscViewerDestroy(&linesearch->monitor);CHKERRQ(ierr);
566bf7f4e0aSPeter Brune   }
567bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
568bf7f4e0aSPeter Brune }
569bf7f4e0aSPeter Brune 
570bf7f4e0aSPeter Brune #undef __FUNCT__
571f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetMonitor"
572f40b411bSPeter Brune /*@
573*cd7522eaSPeter Brune    SNESLineSearchGetMonitor - Gets the PetscViewer instance for the line search monitor.
5746a388c36SPeter Brune 
575f40b411bSPeter Brune    Input Parameters:
576f40b411bSPeter Brune .  linesearch - linesearch context.
577f40b411bSPeter Brune 
578f40b411bSPeter Brune    Input Parameters:
579f40b411bSPeter Brune .  monitor - monitor context.
580f40b411bSPeter Brune 
581f40b411bSPeter Brune    Logically Collective on SNES
582f40b411bSPeter Brune 
583f40b411bSPeter Brune 
584f40b411bSPeter Brune    Options Database Keys:
585*cd7522eaSPeter Brune .   -snes_linesearch_monitor - enables the monitor.
586f40b411bSPeter Brune 
587f40b411bSPeter Brune    Level: intermediate
588f40b411bSPeter Brune 
589f40b411bSPeter Brune 
590*cd7522eaSPeter Brune .seealso: SNESLineSearchSetMonitor(), PetscViewer
591f40b411bSPeter Brune @*/
592f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetMonitor(SNESLineSearch linesearch, PetscViewer *monitor)
5936a388c36SPeter Brune {
5946a388c36SPeter Brune 
5956a388c36SPeter Brune   PetscFunctionBegin;
596f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
5976a388c36SPeter Brune   if (monitor) {
5986a388c36SPeter Brune     PetscValidPointer(monitor, 2);
5996a388c36SPeter Brune     *monitor = linesearch->monitor;
6006a388c36SPeter Brune   }
6016a388c36SPeter Brune   PetscFunctionReturn(0);
6026a388c36SPeter Brune }
6036a388c36SPeter Brune 
6046a388c36SPeter Brune #undef __FUNCT__
605f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetFromOptions"
606f40b411bSPeter Brune /*@
607f1c6b773SPeter Brune    SNESLineSearchSetFromOptions - Sets options for the line search
608f40b411bSPeter Brune 
609f40b411bSPeter Brune    Input Parameters:
610f40b411bSPeter Brune .  linesearch - linesearch context.
611f40b411bSPeter Brune 
612f40b411bSPeter Brune    Options Database Keys:
613*cd7522eaSPeter Brune + -snes_linesearch_type - The Line search method
614*cd7522eaSPeter Brune . -snes_linesearch_monitor - Print progress of line searches
615*cd7522eaSPeter Brune . -snes_linesearch_damping - The linesearch damping parameter.
616*cd7522eaSPeter Brune . -snes_linesearch_norms   - Turn on/off the linesearch norms
617*cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess.
618*cd7522eaSPeter Brune - -snes_linesearch_max_it - The number of iterations for iterative line searches.
619f40b411bSPeter Brune 
620f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
621f40b411bSPeter Brune 
622f40b411bSPeter Brune    Level: intermediate
623f40b411bSPeter Brune 
624f40b411bSPeter Brune 
625f1c6b773SPeter Brune .seealso: SNESLineSearchCreate()
626f40b411bSPeter Brune @*/
627f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch) {
628bf7f4e0aSPeter Brune   PetscErrorCode ierr;
629f1c6b773SPeter Brune   const char     *deft = SNES_LINESEARCH_BASIC;
630bf7f4e0aSPeter Brune   char           type[256];
631bf7f4e0aSPeter Brune   PetscBool      flg, set;
632bf7f4e0aSPeter Brune   PetscFunctionBegin;
633f1c6b773SPeter Brune   if (!SNESLineSearchRegisterAllCalled) {ierr = SNESLineSearchRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
634bf7f4e0aSPeter Brune 
635bf7f4e0aSPeter Brune   ierr = PetscObjectOptionsBegin((PetscObject)linesearch);CHKERRQ(ierr);
636bf7f4e0aSPeter Brune   if (((PetscObject)linesearch)->type_name) {
637bf7f4e0aSPeter Brune     deft = ((PetscObject)linesearch)->type_name;
638bf7f4e0aSPeter Brune   }
6397a35526eSPeter Brune   ierr = PetscOptionsList("-snes_linesearch_type","Line-search method","SNESLineSearchSetType",SNESLineSearchList,deft,type,256,&flg);CHKERRQ(ierr);
640bf7f4e0aSPeter Brune   if (flg) {
641f1c6b773SPeter Brune     ierr = SNESLineSearchSetType(linesearch,type);CHKERRQ(ierr);
642bf7f4e0aSPeter Brune   } else if (!((PetscObject)linesearch)->type_name) {
643f1c6b773SPeter Brune     ierr = SNESLineSearchSetType(linesearch,deft);CHKERRQ(ierr);
644bf7f4e0aSPeter Brune   }
645bf7f4e0aSPeter Brune   if (linesearch->ops->setfromoptions) {
646bf7f4e0aSPeter Brune     (*linesearch->ops->setfromoptions)(linesearch);CHKERRQ(ierr);
647bf7f4e0aSPeter Brune   }
648bf7f4e0aSPeter Brune 
6497a35526eSPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_monitor","Print progress of line searches","SNESSNESLineSearchSetMonitor",
650bf7f4e0aSPeter Brune                           linesearch->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
651f1c6b773SPeter Brune   if (set) {ierr = SNESLineSearchSetMonitor(linesearch,flg);CHKERRQ(ierr);}
652bf7f4e0aSPeter Brune 
6537a35526eSPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_damping","Line search damping and initial step guess","SNESLineSearchSetDamping",linesearch->damping,&linesearch->damping,0);CHKERRQ(ierr);
6547a35526eSPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_rtol","Tolerance for iterative line search","SNESLineSearchSetRTolerance",linesearch->rtol,&linesearch->rtol,0);CHKERRQ(ierr);
6557a35526eSPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_norms","Compute final norms in line search","SNESLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,0);CHKERRQ(ierr);
6567a35526eSPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_keeplambda","Use previous lambda as damping","SNESLineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,0);CHKERRQ(ierr);
6577a35526eSPeter Brune   ierr = PetscOptionsInt("-snes_linesearch_max_it","Maximum iterations for iterative line searches","",linesearch->max_its,&linesearch->max_its,0);CHKERRQ(ierr);
6587a35526eSPeter Brune 
6597a35526eSPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_precheck_picard","Use a correction that sometimes improves convergence of Picard iteration","SNESLineSearchPreCheckPicard",flg,&flg,&set);CHKERRQ(ierr);
6607a35526eSPeter Brune   if (set) {
6617a35526eSPeter Brune     if (flg) {
6627a35526eSPeter Brune       linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */
6637a35526eSPeter Brune       ierr = PetscOptionsReal("-snes_linesearch_precheck_picard_angle","Maximum angle at which to activate the correction",
6647a35526eSPeter Brune                               "none",linesearch->precheck_picard_angle,&linesearch->precheck_picard_angle,PETSC_NULL);CHKERRQ(ierr);
6657a35526eSPeter Brune       ierr = SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&linesearch->precheck_picard_angle);CHKERRQ(ierr);
6667a35526eSPeter Brune     } else {
6677a35526eSPeter Brune       ierr = SNESLineSearchSetPreCheck(linesearch,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
6687a35526eSPeter Brune     }
6697a35526eSPeter Brune   }
6707a35526eSPeter Brune 
671bf7f4e0aSPeter Brune   ierr = PetscObjectProcessOptionsHandlers((PetscObject)linesearch);CHKERRQ(ierr);
672bf7f4e0aSPeter Brune   ierr = PetscOptionsEnd();CHKERRQ(ierr);
673bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
674bf7f4e0aSPeter Brune }
675bf7f4e0aSPeter Brune 
676bf7f4e0aSPeter Brune #undef __FUNCT__
677f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchView"
678f40b411bSPeter Brune /*@
679*cd7522eaSPeter Brune    SNESLineSearchView - Prints useful information about the line search not
680*cd7522eaSPeter Brune    related to an individual call.
681f40b411bSPeter Brune 
682f40b411bSPeter Brune    Input Parameters:
683f40b411bSPeter Brune .  linesearch - linesearch context.
684f40b411bSPeter Brune 
685f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
686f40b411bSPeter Brune 
687f40b411bSPeter Brune    Level: intermediate
688f40b411bSPeter Brune 
689f1c6b773SPeter Brune .seealso: SNESLineSearchCreate()
690f40b411bSPeter Brune @*/
691f1c6b773SPeter Brune PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch) {
692bf7f4e0aSPeter Brune   PetscFunctionBegin;
693f40b411bSPeter Brune 
694bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
695bf7f4e0aSPeter Brune }
696bf7f4e0aSPeter Brune 
697bf7f4e0aSPeter Brune #undef __FUNCT__
698f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetType"
699ea5d4fccSPeter Brune /*@C
700f1c6b773SPeter Brune    SNESLineSearchSetType - Sets the linesearch type
701f40b411bSPeter Brune 
702f40b411bSPeter Brune    Input Parameters:
703f40b411bSPeter Brune +  linesearch - linesearch context.
704f40b411bSPeter Brune -  type - The type of line search to be used
705f40b411bSPeter Brune 
706*cd7522eaSPeter Brune    Available Types:
707*cd7522eaSPeter Brune +  basic - Simple damping line search.
708*cd7522eaSPeter Brune .  bt - Backtracking line search over the L2 norm of the function.
709*cd7522eaSPeter Brune .  l2 - Secant line search over the L2 norm of the function.
710*cd7522eaSPeter Brune .  cp - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x).
711*cd7522eaSPeter Brune -  shell - User provided SNESLineSearch implementation.
712*cd7522eaSPeter Brune 
713f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
714f40b411bSPeter Brune 
715f40b411bSPeter Brune    Level: intermediate
716f40b411bSPeter Brune 
717f40b411bSPeter Brune 
718f1c6b773SPeter Brune .seealso: SNESLineSearchCreate()
719f40b411bSPeter Brune @*/
720f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, const SNESLineSearchType type)
721bf7f4e0aSPeter Brune {
722bf7f4e0aSPeter Brune 
723f1c6b773SPeter Brune   PetscErrorCode ierr,(*r)(SNESLineSearch);
724bf7f4e0aSPeter Brune   PetscBool      match;
725bf7f4e0aSPeter Brune 
726bf7f4e0aSPeter Brune   PetscFunctionBegin;
727f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
728bf7f4e0aSPeter Brune   PetscValidCharPointer(type,2);
729bf7f4e0aSPeter Brune 
730bf7f4e0aSPeter Brune   ierr = PetscTypeCompare((PetscObject)linesearch,type,&match);CHKERRQ(ierr);
731bf7f4e0aSPeter Brune   if (match) PetscFunctionReturn(0);
732bf7f4e0aSPeter Brune 
733f1c6b773SPeter Brune   ierr =  PetscFListFind(SNESLineSearchList,((PetscObject)linesearch)->comm,type,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr);
734bf7f4e0aSPeter Brune   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type);
735bf7f4e0aSPeter Brune   /* Destroy the previous private linesearch context */
736bf7f4e0aSPeter Brune   if (linesearch->ops->destroy) {
737bf7f4e0aSPeter Brune     ierr = (*(linesearch)->ops->destroy)(linesearch);CHKERRQ(ierr);
738bf7f4e0aSPeter Brune     linesearch->ops->destroy = PETSC_NULL;
739bf7f4e0aSPeter Brune   }
740f1c6b773SPeter Brune   /* Reinitialize function pointers in SNESLineSearchOps structure */
741bf7f4e0aSPeter Brune   linesearch->ops->apply          = 0;
742bf7f4e0aSPeter Brune   linesearch->ops->view           = 0;
743bf7f4e0aSPeter Brune   linesearch->ops->setfromoptions = 0;
744bf7f4e0aSPeter Brune   linesearch->ops->destroy        = 0;
745bf7f4e0aSPeter Brune 
746bf7f4e0aSPeter Brune   ierr = PetscObjectChangeTypeName((PetscObject)linesearch,type);CHKERRQ(ierr);
747bf7f4e0aSPeter Brune   ierr = (*r)(linesearch);CHKERRQ(ierr);
748bf7f4e0aSPeter Brune #if defined(PETSC_HAVE_AMS)
749bf7f4e0aSPeter Brune   if (PetscAMSPublishAll) {
750bf7f4e0aSPeter Brune     ierr = PetscObjectAMSPublish((PetscObject)linesearch);CHKERRQ(ierr);
751bf7f4e0aSPeter Brune   }
752bf7f4e0aSPeter Brune #endif
753bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
754bf7f4e0aSPeter Brune }
755bf7f4e0aSPeter Brune 
756bf7f4e0aSPeter Brune #undef __FUNCT__
757f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetSNES"
758f40b411bSPeter Brune /*@
759f1c6b773SPeter Brune    SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation
760f40b411bSPeter Brune 
761f40b411bSPeter Brune    Input Parameters:
762f40b411bSPeter Brune +  linesearch - linesearch context.
763f40b411bSPeter Brune -  snes - The snes instance
764f40b411bSPeter Brune 
765f40b411bSPeter Brune    Level: intermediate
766f40b411bSPeter Brune 
767f40b411bSPeter Brune 
768*cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
769f40b411bSPeter Brune @*/
770f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes){
771bf7f4e0aSPeter Brune   PetscFunctionBegin;
772f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
773bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,2);
774bf7f4e0aSPeter Brune   linesearch->snes = snes;
775bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
776bf7f4e0aSPeter Brune }
777bf7f4e0aSPeter Brune 
778bf7f4e0aSPeter Brune #undef __FUNCT__
779f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetSNES"
780f40b411bSPeter Brune /*@
781f1c6b773SPeter Brune    SNESLineSearchGetSNES - Gets the SNES for the linesearch for function evaluation
782f40b411bSPeter Brune 
783f40b411bSPeter Brune    Input Parameters:
784f40b411bSPeter Brune .  linesearch - linesearch context.
785f40b411bSPeter Brune 
786f40b411bSPeter Brune    Output Parameters:
787f40b411bSPeter Brune .  snes - The snes instance
788f40b411bSPeter Brune 
789f40b411bSPeter Brune    Level: intermediate
790f40b411bSPeter Brune 
791*cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
792f40b411bSPeter Brune @*/
793f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes){
794bf7f4e0aSPeter Brune   PetscFunctionBegin;
795f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
7966a388c36SPeter Brune   PetscValidPointer(snes, 2);
797bf7f4e0aSPeter Brune   *snes = linesearch->snes;
798bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
799bf7f4e0aSPeter Brune }
800bf7f4e0aSPeter Brune 
8016a388c36SPeter Brune #undef __FUNCT__
802f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetLambda"
803f40b411bSPeter Brune /*@
804f1c6b773SPeter Brune    SNESLineSearchGetLambda - Gets the last linesearch steplength discovered.
805f40b411bSPeter Brune 
806f40b411bSPeter Brune    Input Parameters:
807f40b411bSPeter Brune .  linesearch - linesearch context.
808f40b411bSPeter Brune 
809f40b411bSPeter Brune    Output Parameters:
810*cd7522eaSPeter Brune .  lambda - The last steplength computed during SNESLineSearchApply()
811f40b411bSPeter Brune 
812f40b411bSPeter Brune    Level: intermediate
813f40b411bSPeter Brune 
814*cd7522eaSPeter Brune .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetDamping(), SNESLineSearchApply()
815f40b411bSPeter Brune @*/
816f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda)
8176a388c36SPeter Brune {
8186a388c36SPeter Brune   PetscFunctionBegin;
819f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
8206a388c36SPeter Brune   PetscValidPointer(lambda, 2);
8216a388c36SPeter Brune   *lambda = linesearch->lambda;
8226a388c36SPeter Brune   PetscFunctionReturn(0);
8236a388c36SPeter Brune }
8246a388c36SPeter Brune 
8256a388c36SPeter Brune #undef __FUNCT__
826f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetLambda"
827f40b411bSPeter Brune /*@
828f1c6b773SPeter Brune    SNESLineSearchSetLambda - Sets the linesearch steplength.
829f40b411bSPeter Brune 
830f40b411bSPeter Brune    Input Parameters:
831f40b411bSPeter Brune +  linesearch - linesearch context.
832f40b411bSPeter Brune -  lambda - The last steplength.
833f40b411bSPeter Brune 
834*cd7522eaSPeter Brune    Notes:
835*cd7522eaSPeter Brune    This routine is typically used within implementations of SNESLineSearchApply
836*cd7522eaSPeter Brune    to set the final steplength.  This routine (and SNESLineSearchGetLambda()) were
837*cd7522eaSPeter Brune    added in order to facilitate Quasi-Newton methods that use the previous steplength
838*cd7522eaSPeter Brune    as an inner scaling parameter.
839*cd7522eaSPeter Brune 
840f40b411bSPeter Brune    Level: intermediate
841f40b411bSPeter Brune 
842f1c6b773SPeter Brune .seealso: SNESLineSearchGetLambda()
843f40b411bSPeter Brune @*/
844f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
8456a388c36SPeter Brune {
8466a388c36SPeter Brune   PetscFunctionBegin;
847f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
8486a388c36SPeter Brune   linesearch->lambda = lambda;
8496a388c36SPeter Brune   PetscFunctionReturn(0);
8506a388c36SPeter Brune }
8516a388c36SPeter Brune 
8526a388c36SPeter Brune #undef  __FUNCT__
853f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetTolerances"
854f40b411bSPeter Brune /*@
855f1c6b773SPeter Brune    SNESLineSearchGetTolerances - Gets the tolerances for the method
856f40b411bSPeter Brune 
857f40b411bSPeter Brune    Input Parameters:
858f40b411bSPeter Brune .  linesearch - linesearch context.
859f40b411bSPeter Brune 
860f40b411bSPeter Brune    Output Parameters:
861516fe3c3SPeter Brune +  steptol - The minimum steplength
862516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
863516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
864516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
865516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
866f40b411bSPeter Brune 
867516fe3c3SPeter Brune    Level: advanced
868516fe3c3SPeter Brune 
869f1c6b773SPeter Brune .seealso: SNESLineSearchSetTolerances()
870f40b411bSPeter Brune @*/
871f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
8726a388c36SPeter Brune {
8736a388c36SPeter Brune   PetscFunctionBegin;
874f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
875516fe3c3SPeter Brune   if (steptol) {
8766a388c36SPeter Brune     PetscValidPointer(steptol, 2);
8776a388c36SPeter Brune     *steptol = linesearch->steptol;
878516fe3c3SPeter Brune   }
879516fe3c3SPeter Brune   if (maxstep) {
880516fe3c3SPeter Brune     PetscValidPointer(maxstep, 3);
881516fe3c3SPeter Brune     *maxstep = linesearch->maxstep;
882516fe3c3SPeter Brune   }
883516fe3c3SPeter Brune   if (rtol) {
884516fe3c3SPeter Brune     PetscValidPointer(rtol, 4);
885516fe3c3SPeter Brune     *rtol = linesearch->rtol;
886516fe3c3SPeter Brune   }
887516fe3c3SPeter Brune   if (atol) {
888516fe3c3SPeter Brune     PetscValidPointer(atol, 5);
889516fe3c3SPeter Brune     *atol = linesearch->atol;
890516fe3c3SPeter Brune   }
891516fe3c3SPeter Brune   if (ltol) {
892516fe3c3SPeter Brune     PetscValidPointer(ltol, 6);
893516fe3c3SPeter Brune     *ltol = linesearch->ltol;
894516fe3c3SPeter Brune   }
895516fe3c3SPeter Brune   if (max_its) {
896516fe3c3SPeter Brune     PetscValidPointer(max_its, 7);
897516fe3c3SPeter Brune     *max_its = linesearch->max_its;
898516fe3c3SPeter Brune   }
8996a388c36SPeter Brune   PetscFunctionReturn(0);
9006a388c36SPeter Brune }
9016a388c36SPeter Brune 
9026a388c36SPeter Brune #undef  __FUNCT__
903f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetTolerances"
904f40b411bSPeter Brune /*@
905f1c6b773SPeter Brune    SNESLineSearchSetTolerances - Sets the tolerances for the method
906f40b411bSPeter Brune 
907f40b411bSPeter Brune    Input Parameters:
908516fe3c3SPeter Brune +  linesearch - linesearch context.
909516fe3c3SPeter Brune .  steptol - The minimum steplength
910516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
911516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
912516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
913516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
914f40b411bSPeter Brune 
915f40b411bSPeter Brune 
916516fe3c3SPeter Brune    Level: advanced
917516fe3c3SPeter Brune 
918f1c6b773SPeter Brune .seealso: SNESLineSearchGetTolerances()
919f40b411bSPeter Brune @*/
920f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its)
9216a388c36SPeter Brune {
9226a388c36SPeter Brune   PetscFunctionBegin;
923f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
9246a388c36SPeter Brune   linesearch->steptol = steptol;
925516fe3c3SPeter Brune   linesearch->maxstep = maxstep;
926516fe3c3SPeter Brune   linesearch->rtol = rtol;
927516fe3c3SPeter Brune   linesearch->atol = atol;
928516fe3c3SPeter Brune   linesearch->ltol = ltol;
929516fe3c3SPeter Brune   linesearch->max_its = max_its;
9306a388c36SPeter Brune   PetscFunctionReturn(0);
9316a388c36SPeter Brune }
9326a388c36SPeter Brune 
933516fe3c3SPeter Brune 
9346a388c36SPeter Brune #undef __FUNCT__
935f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetDamping"
936f40b411bSPeter Brune /*@
937f1c6b773SPeter Brune    SNESLineSearchGetDamping - Gets the line search damping parameter.
938f40b411bSPeter Brune 
939f40b411bSPeter Brune    Input Parameters:
940f40b411bSPeter Brune .  linesearch - linesearch context.
941f40b411bSPeter Brune 
942f40b411bSPeter Brune    Output Parameters:
943f40b411bSPeter Brune .  damping - The damping parameter.
944f40b411bSPeter Brune 
945f40b411bSPeter Brune    Level: intermediate
946f40b411bSPeter Brune 
947f1c6b773SPeter Brune .seealso: SNESLineSearchGetStepTolerance()
948f40b411bSPeter Brune @*/
949f40b411bSPeter Brune 
950f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping)
9516a388c36SPeter Brune {
9526a388c36SPeter Brune   PetscFunctionBegin;
953f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
9546a388c36SPeter Brune   PetscValidPointer(damping, 2);
9556a388c36SPeter Brune   *damping = linesearch->damping;
9566a388c36SPeter Brune   PetscFunctionReturn(0);
9576a388c36SPeter Brune }
9586a388c36SPeter Brune 
9596a388c36SPeter Brune #undef __FUNCT__
960f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetDamping"
961f40b411bSPeter Brune /*@
962f1c6b773SPeter Brune    SNESLineSearchSetDamping - Sets the line search damping paramter.
963f40b411bSPeter Brune 
964f40b411bSPeter Brune    Input Parameters:
965f40b411bSPeter Brune .  linesearch - linesearch context.
966f40b411bSPeter Brune .  damping - The damping parameter.
967f40b411bSPeter Brune 
968f40b411bSPeter Brune    Level: intermediate
969f40b411bSPeter Brune 
970*cd7522eaSPeter Brune    Notes:
971*cd7522eaSPeter Brune    The basic line search merely takes the update step scaled by the damping parameter.
972*cd7522eaSPeter Brune    The use of the damping parameter in the l2 and cp line searches is much more subtle;
973*cd7522eaSPeter Brune    it is used as a starting point in calculating the secant step; however, the eventual
974*cd7522eaSPeter Brune    step may be of greater length than the damping parameter.  In the bt line search it is
975*cd7522eaSPeter Brune    used as the maximum possible step length, as the bt line search only backtracks.
976*cd7522eaSPeter Brune 
977f1c6b773SPeter Brune .seealso: SNESLineSearchGetDamping()
978f40b411bSPeter Brune @*/
979f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping)
9806a388c36SPeter Brune {
9816a388c36SPeter Brune   PetscFunctionBegin;
982f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
9836a388c36SPeter Brune   linesearch->damping = damping;
9846a388c36SPeter Brune   PetscFunctionReturn(0);
9856a388c36SPeter Brune }
9866a388c36SPeter Brune 
9876a388c36SPeter Brune #undef __FUNCT__
988f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetNorms"
989f40b411bSPeter Brune /*@
990f1c6b773SPeter Brune    SNESLineSearchGetNorms - Gets the norms for for X, Y, and F.
991f40b411bSPeter Brune 
992f40b411bSPeter Brune    Input Parameters:
993f40b411bSPeter Brune .  linesearch - linesearch context.
994f40b411bSPeter Brune 
995f40b411bSPeter Brune    Output Parameters:
996f40b411bSPeter Brune +  xnorm - The norm of the current solution
997f40b411bSPeter Brune .  fnorm - The norm of the current function
998f40b411bSPeter Brune -  ynorm - The norm of the current update
999f40b411bSPeter Brune 
1000*cd7522eaSPeter Brune    Notes:
1001*cd7522eaSPeter Brune    This function is mainly called from SNES implementations.
1002*cd7522eaSPeter Brune 
1003f40b411bSPeter Brune    Level: intermediate
1004f40b411bSPeter Brune 
1005f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms() SNESLineSearchGetVecs()
1006f40b411bSPeter Brune @*/
1007f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm)
1008bf7f4e0aSPeter Brune {
1009bf7f4e0aSPeter Brune   PetscFunctionBegin;
1010f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1011bf7f4e0aSPeter Brune   if (xnorm) {
1012bf7f4e0aSPeter Brune     *xnorm = linesearch->xnorm;
1013bf7f4e0aSPeter Brune   }
1014bf7f4e0aSPeter Brune   if (fnorm) {
1015bf7f4e0aSPeter Brune     *fnorm = linesearch->fnorm;
1016bf7f4e0aSPeter Brune   }
1017bf7f4e0aSPeter Brune   if (ynorm) {
1018bf7f4e0aSPeter Brune     *ynorm = linesearch->ynorm;
1019bf7f4e0aSPeter Brune   }
1020bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1021bf7f4e0aSPeter Brune }
1022bf7f4e0aSPeter Brune 
1023e7058c64SPeter Brune #undef __FUNCT__
1024f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetNorms"
1025f40b411bSPeter Brune /*@
1026f1c6b773SPeter Brune    SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F.
1027f40b411bSPeter Brune 
1028f40b411bSPeter Brune    Input Parameters:
1029f40b411bSPeter Brune +  linesearch - linesearch context.
1030f40b411bSPeter Brune .  xnorm - The norm of the current solution
1031f40b411bSPeter Brune .  fnorm - The norm of the current function
1032f40b411bSPeter Brune -  ynorm - The norm of the current update
1033f40b411bSPeter Brune 
1034f40b411bSPeter Brune    Level: intermediate
1035f40b411bSPeter Brune 
1036f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1037f40b411bSPeter Brune @*/
1038f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
10396a388c36SPeter Brune {
10406a388c36SPeter Brune   PetscFunctionBegin;
1041f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
10426a388c36SPeter Brune   linesearch->xnorm = xnorm;
10436a388c36SPeter Brune   linesearch->fnorm = fnorm;
10446a388c36SPeter Brune   linesearch->ynorm = ynorm;
10456a388c36SPeter Brune   PetscFunctionReturn(0);
10466a388c36SPeter Brune }
10476a388c36SPeter Brune 
10486a388c36SPeter Brune #undef __FUNCT__
1049f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchComputeNorms"
1050f40b411bSPeter Brune /*@
1051f1c6b773SPeter Brune    SNESLineSearchComputeNorms - Computes the norms of X, F, and Y.
1052f40b411bSPeter Brune 
1053f40b411bSPeter Brune    Input Parameters:
1054f40b411bSPeter Brune .  linesearch - linesearch context.
1055f40b411bSPeter Brune 
1056f40b411bSPeter Brune    Options Database Keys:
1057*cd7522eaSPeter Brune .   -snes_linesearch_norms - turn norm computation on or off.
1058f40b411bSPeter Brune 
1059f40b411bSPeter Brune    Level: intermediate
1060f40b411bSPeter Brune 
1061f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms, SNESLineSearchSetNorms()
1062f40b411bSPeter Brune @*/
1063f1c6b773SPeter Brune PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
10646a388c36SPeter Brune {
10656a388c36SPeter Brune   PetscErrorCode ierr;
10669bd66eb0SPeter Brune   SNES snes;
10676a388c36SPeter Brune   PetscFunctionBegin;
10686a388c36SPeter Brune   if (linesearch->norms) {
10699bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
1070f1c6b773SPeter Brune       ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
10719bd66eb0SPeter Brune       ierr = VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
10729bd66eb0SPeter Brune       ierr = VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
10739bd66eb0SPeter Brune       ierr = (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);CHKERRQ(ierr);
10749bd66eb0SPeter Brune     } else {
10756a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_func,   NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
10766a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_sol,    NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
10776a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
10786a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_func,     NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
10796a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_sol,      NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
10806a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_update,   NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
10816a388c36SPeter Brune     }
10829bd66eb0SPeter Brune   }
10836a388c36SPeter Brune   PetscFunctionReturn(0);
10846a388c36SPeter Brune }
10856a388c36SPeter Brune 
10866a388c36SPeter Brune #undef __FUNCT__
1087f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetVecs"
1088f40b411bSPeter Brune /*@
1089f1c6b773SPeter Brune    SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context
1090f40b411bSPeter Brune 
1091f40b411bSPeter Brune    Input Parameters:
1092f40b411bSPeter Brune .  linesearch - linesearch context.
1093f40b411bSPeter Brune 
1094f40b411bSPeter Brune    Output Parameters:
1095f40b411bSPeter Brune +  X - The old solution
1096f40b411bSPeter Brune .  F - The old function
1097f40b411bSPeter Brune .  Y - The search direction
1098f40b411bSPeter Brune .  W - The new solution
1099f40b411bSPeter Brune -  G - The new function
1100f40b411bSPeter Brune 
1101f40b411bSPeter Brune    Level: intermediate
1102f40b411bSPeter Brune 
1103f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1104f40b411bSPeter Brune @*/
1105f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G) {
11066a388c36SPeter Brune   PetscFunctionBegin;
1107f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
11086a388c36SPeter Brune   if (X) {
11096a388c36SPeter Brune     PetscValidPointer(X, 2);
11106a388c36SPeter Brune     *X = linesearch->vec_sol;
11116a388c36SPeter Brune   }
11126a388c36SPeter Brune   if (F) {
11136a388c36SPeter Brune     PetscValidPointer(F, 3);
11146a388c36SPeter Brune     *F = linesearch->vec_func;
11156a388c36SPeter Brune   }
11166a388c36SPeter Brune   if (Y) {
11176a388c36SPeter Brune     PetscValidPointer(Y, 4);
11186a388c36SPeter Brune     *Y = linesearch->vec_update;
11196a388c36SPeter Brune   }
11206a388c36SPeter Brune   if (W) {
11216a388c36SPeter Brune     PetscValidPointer(W, 5);
11226a388c36SPeter Brune     *W = linesearch->vec_sol_new;
11236a388c36SPeter Brune   }
11246a388c36SPeter Brune   if (G) {
11256a388c36SPeter Brune     PetscValidPointer(G, 6);
11266a388c36SPeter Brune     *G = linesearch->vec_func_new;
11276a388c36SPeter Brune   }
11286a388c36SPeter Brune 
11296a388c36SPeter Brune   PetscFunctionReturn(0);
11306a388c36SPeter Brune }
11316a388c36SPeter Brune 
11326a388c36SPeter Brune #undef __FUNCT__
1133f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetVecs"
1134f40b411bSPeter Brune /*@
1135f1c6b773SPeter Brune    SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context
1136f40b411bSPeter Brune 
1137f40b411bSPeter Brune    Input Parameters:
1138f40b411bSPeter Brune +  linesearch - linesearch context.
1139f40b411bSPeter Brune .  X - The old solution
1140f40b411bSPeter Brune .  F - The old function
1141f40b411bSPeter Brune .  Y - The search direction
1142f40b411bSPeter Brune .  W - The new solution
1143f40b411bSPeter Brune -  G - The new function
1144f40b411bSPeter Brune 
1145f40b411bSPeter Brune    Level: intermediate
1146f40b411bSPeter Brune 
1147f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms(), SNESLineSearchGetVecs()
1148f40b411bSPeter Brune @*/
1149f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G) {
11506a388c36SPeter Brune   PetscFunctionBegin;
1151f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
11526a388c36SPeter Brune   if (X) {
11536a388c36SPeter Brune     PetscValidHeaderSpecific(X,VEC_CLASSID,2);
11546a388c36SPeter Brune     linesearch->vec_sol = X;
11556a388c36SPeter Brune   }
11566a388c36SPeter Brune   if (F) {
11576a388c36SPeter Brune     PetscValidHeaderSpecific(F,VEC_CLASSID,3);
11586a388c36SPeter Brune     linesearch->vec_func = F;
11596a388c36SPeter Brune   }
11606a388c36SPeter Brune   if (Y) {
11616a388c36SPeter Brune     PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
11626a388c36SPeter Brune     linesearch->vec_update = Y;
11636a388c36SPeter Brune   }
11646a388c36SPeter Brune   if (W) {
11656a388c36SPeter Brune     PetscValidHeaderSpecific(W,VEC_CLASSID,5);
11666a388c36SPeter Brune     linesearch->vec_sol_new = W;
11676a388c36SPeter Brune   }
11686a388c36SPeter Brune   if (G) {
11696a388c36SPeter Brune     PetscValidHeaderSpecific(G,VEC_CLASSID,6);
11706a388c36SPeter Brune     linesearch->vec_func_new = G;
11716a388c36SPeter Brune   }
11726a388c36SPeter Brune 
11736a388c36SPeter Brune   PetscFunctionReturn(0);
11746a388c36SPeter Brune }
11756a388c36SPeter Brune 
11766a388c36SPeter Brune #undef __FUNCT__
1177f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchAppendOptionsPrefix"
1178e7058c64SPeter Brune /*@C
1179f1c6b773SPeter Brune    SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1180e7058c64SPeter Brune    SNES options in the database.
1181e7058c64SPeter Brune 
1182*cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
1183e7058c64SPeter Brune 
1184e7058c64SPeter Brune    Input Parameters:
1185e7058c64SPeter Brune +  snes - the SNES context
1186e7058c64SPeter Brune -  prefix - the prefix to prepend to all option names
1187e7058c64SPeter Brune 
1188e7058c64SPeter Brune    Notes:
1189e7058c64SPeter Brune    A hyphen (-) must NOT be given at the beginning of the prefix name.
1190e7058c64SPeter Brune    The first character of all runtime options is AUTOMATICALLY the hyphen.
1191e7058c64SPeter Brune 
1192e7058c64SPeter Brune    Level: advanced
1193e7058c64SPeter Brune 
1194f1c6b773SPeter Brune .keywords: SNESLineSearch, append, options, prefix, database
1195e7058c64SPeter Brune 
1196e7058c64SPeter Brune .seealso: SNESGetOptionsPrefix()
1197e7058c64SPeter Brune @*/
1198f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[])
1199e7058c64SPeter Brune {
1200e7058c64SPeter Brune   PetscErrorCode ierr;
1201e7058c64SPeter Brune 
1202e7058c64SPeter Brune   PetscFunctionBegin;
1203f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1204e7058c64SPeter Brune   ierr = PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1205e7058c64SPeter Brune   PetscFunctionReturn(0);
1206e7058c64SPeter Brune }
1207e7058c64SPeter Brune 
1208e7058c64SPeter Brune #undef __FUNCT__
1209f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetOptionsPrefix"
1210e7058c64SPeter Brune /*@C
1211f1c6b773SPeter Brune    SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all
1212f1c6b773SPeter Brune    SNESLineSearch options in the database.
1213e7058c64SPeter Brune 
1214e7058c64SPeter Brune    Not Collective
1215e7058c64SPeter Brune 
1216e7058c64SPeter Brune    Input Parameter:
1217*cd7522eaSPeter Brune .  linesearch - the SNESLineSearch context
1218e7058c64SPeter Brune 
1219e7058c64SPeter Brune    Output Parameter:
1220e7058c64SPeter Brune .  prefix - pointer to the prefix string used
1221e7058c64SPeter Brune 
1222e7058c64SPeter Brune    Notes: On the fortran side, the user should pass in a string 'prefix' of
1223e7058c64SPeter Brune    sufficient length to hold the prefix.
1224e7058c64SPeter Brune 
1225e7058c64SPeter Brune    Level: advanced
1226e7058c64SPeter Brune 
1227f1c6b773SPeter Brune .keywords: SNESLineSearch, get, options, prefix, database
1228e7058c64SPeter Brune 
1229e7058c64SPeter Brune .seealso: SNESAppendOptionsPrefix()
1230e7058c64SPeter Brune @*/
1231f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[])
1232e7058c64SPeter Brune {
1233e7058c64SPeter Brune   PetscErrorCode ierr;
1234e7058c64SPeter Brune 
1235e7058c64SPeter Brune   PetscFunctionBegin;
1236f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1237e7058c64SPeter Brune   ierr = PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1238e7058c64SPeter Brune   PetscFunctionReturn(0);
1239e7058c64SPeter Brune }
1240bf7f4e0aSPeter Brune 
1241bf7f4e0aSPeter Brune #undef __FUNCT__
1242f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetWork"
1243f40b411bSPeter Brune /*@
1244f1c6b773SPeter Brune    SNESLineSearchGetWork - Gets work vectors for the line search.
1245f40b411bSPeter Brune 
1246f40b411bSPeter Brune    Input Parameter:
1247f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
1248f40b411bSPeter Brune -  nwork - the number of work vectors
1249f40b411bSPeter Brune 
1250f40b411bSPeter Brune    Level: developer
1251f40b411bSPeter Brune 
1252*cd7522eaSPeter Brune    Notes:
1253*cd7522eaSPeter Brune    This is typically called at the beginning of a SNESLineSearch or SNESLineSearchShell implementation.
1254*cd7522eaSPeter Brune 
1255f1c6b773SPeter Brune .keywords: SNESLineSearch, work, vector
1256f40b411bSPeter Brune 
1257f40b411bSPeter Brune .seealso: SNESDefaultGetWork()
1258f40b411bSPeter Brune @*/
1259f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetWork(SNESLineSearch linesearch, PetscInt nwork)
1260bf7f4e0aSPeter Brune {
1261bf7f4e0aSPeter Brune   PetscErrorCode ierr;
1262bf7f4e0aSPeter Brune   PetscFunctionBegin;
1263bf7f4e0aSPeter Brune   if (linesearch->vec_sol) {
1264bf7f4e0aSPeter Brune     ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr);
1265bf7f4e0aSPeter Brune   } else {
1266bf7f4e0aSPeter Brune     SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
1267bf7f4e0aSPeter Brune   }
1268bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1269bf7f4e0aSPeter Brune }
1270bf7f4e0aSPeter Brune 
12716a388c36SPeter Brune 
1272bf7f4e0aSPeter Brune #undef __FUNCT__
1273f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetSuccess"
1274f40b411bSPeter Brune /*@
1275f1c6b773SPeter Brune    SNESLineSearchGetSuccess - Gets the success/failure status of the last line search application
1276f40b411bSPeter Brune 
1277f40b411bSPeter Brune    Input Parameters:
1278f40b411bSPeter Brune .  linesearch - linesearch context.
1279f40b411bSPeter Brune 
1280f40b411bSPeter Brune    Output Parameters:
1281f40b411bSPeter Brune .  success - The success or failure status.
1282f40b411bSPeter Brune 
1283*cd7522eaSPeter Brune    Notes:
1284*cd7522eaSPeter Brune    This is typically called after SNESLineSearchApply in order to determine if the line-search failed
1285*cd7522eaSPeter Brune    (and set the SNES convergence accordingly).
1286*cd7522eaSPeter Brune 
1287f40b411bSPeter Brune    Level: intermediate
1288f40b411bSPeter Brune 
1289f1c6b773SPeter Brune .seealso: SNESLineSearchSetSuccess()
1290f40b411bSPeter Brune @*/
1291f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetSuccess(SNESLineSearch linesearch, PetscBool *success)
1292bf7f4e0aSPeter Brune {
1293bf7f4e0aSPeter Brune   PetscFunctionBegin;
1294f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
12956a388c36SPeter Brune   PetscValidPointer(success, 2);
1296bf7f4e0aSPeter Brune   if (success) {
1297bf7f4e0aSPeter Brune     *success = linesearch->success;
1298bf7f4e0aSPeter Brune   }
1299bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1300bf7f4e0aSPeter Brune }
1301bf7f4e0aSPeter Brune 
1302bf7f4e0aSPeter Brune #undef __FUNCT__
1303f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetSuccess"
1304f40b411bSPeter Brune /*@
1305f1c6b773SPeter Brune    SNESLineSearchSetSuccess - Sets the success/failure status of the last line search application
1306f40b411bSPeter Brune 
1307f40b411bSPeter Brune    Input Parameters:
1308f40b411bSPeter Brune +  linesearch - linesearch context.
1309f40b411bSPeter Brune -  success - The success or failure status.
1310f40b411bSPeter Brune 
1311*cd7522eaSPeter Brune    Notes:
1312*cd7522eaSPeter Brune    This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set
1313*cd7522eaSPeter Brune    the success or failure of the line search method.
1314*cd7522eaSPeter Brune 
1315f40b411bSPeter Brune    Level: intermediate
1316f40b411bSPeter Brune 
1317f1c6b773SPeter Brune .seealso: SNESLineSearchGetSuccess()
1318f40b411bSPeter Brune @*/
1319f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetSuccess(SNESLineSearch linesearch, PetscBool success)
13206a388c36SPeter Brune {
1321f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
13226a388c36SPeter Brune   PetscFunctionBegin;
13236a388c36SPeter Brune   linesearch->success = success;
13246a388c36SPeter Brune   PetscFunctionReturn(0);
13256a388c36SPeter Brune }
13266a388c36SPeter Brune 
13276a388c36SPeter Brune #undef __FUNCT__
1328f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetVIFunctions"
13299bd66eb0SPeter Brune /*@C
1330f1c6b773SPeter Brune    SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.
13319bd66eb0SPeter Brune 
13329bd66eb0SPeter Brune    Input Parameters:
13339bd66eb0SPeter Brune +  snes - nonlinear context obtained from SNESCreate()
13349bd66eb0SPeter Brune .  projectfunc - function for projecting the function to the bounds
13359bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
13369bd66eb0SPeter Brune 
13379bd66eb0SPeter Brune    Logically Collective on SNES
13389bd66eb0SPeter Brune 
13399bd66eb0SPeter Brune    Calling sequence of projectfunc:
13409bd66eb0SPeter Brune .vb
13419bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X)
13429bd66eb0SPeter Brune .ve
13439bd66eb0SPeter Brune 
13449bd66eb0SPeter Brune     Input parameters for projectfunc:
13459bd66eb0SPeter Brune +   snes - nonlinear context
13469bd66eb0SPeter Brune -   X - current solution
13479bd66eb0SPeter Brune 
1348*cd7522eaSPeter Brune     Output parameters for projectfunc:
13499bd66eb0SPeter Brune .   X - Projected solution
13509bd66eb0SPeter Brune 
13519bd66eb0SPeter Brune    Calling sequence of normfunc:
13529bd66eb0SPeter Brune .vb
13539bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm)
13549bd66eb0SPeter Brune .ve
13559bd66eb0SPeter Brune 
1356*cd7522eaSPeter Brune     Input parameters for normfunc:
13579bd66eb0SPeter Brune +   snes - nonlinear context
13589bd66eb0SPeter Brune .   X - current solution
13599bd66eb0SPeter Brune -   F - current residual
13609bd66eb0SPeter Brune 
1361*cd7522eaSPeter Brune     Output parameters for normfunc:
13629bd66eb0SPeter Brune .   fnorm - VI-specific norm of the function
13639bd66eb0SPeter Brune 
1364*cd7522eaSPeter Brune     Notes:
1365*cd7522eaSPeter Brune     The VI solvers require projection of the solution to the feasible set.  projectfunc should implement this.
1366*cd7522eaSPeter Brune 
1367*cd7522eaSPeter Brune     The VI solvers require special evaluation of the function norm such that the norm is only calculated
1368*cd7522eaSPeter Brune     on the inactive set.  This should be implemented by normfunc.
13699bd66eb0SPeter Brune 
13709bd66eb0SPeter Brune     Level: developer
13719bd66eb0SPeter Brune 
13729bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, set, line search
13739bd66eb0SPeter Brune 
1374f1c6b773SPeter Brune .seealso: SNESLineSearchGetVIFunctions(), SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck()
13759bd66eb0SPeter Brune @*/
1376f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc)
13779bd66eb0SPeter Brune {
13789bd66eb0SPeter Brune   PetscFunctionBegin;
1379f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
13809bd66eb0SPeter Brune   if (projectfunc) linesearch->ops->viproject = projectfunc;
13819bd66eb0SPeter Brune   if (normfunc) linesearch->ops->vinorm = normfunc;
13829bd66eb0SPeter Brune   PetscFunctionReturn(0);
13839bd66eb0SPeter Brune }
13849bd66eb0SPeter Brune 
13859bd66eb0SPeter Brune #undef __FUNCT__
1386f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetVIFunctions"
13879bd66eb0SPeter Brune /*@C
1388f1c6b773SPeter Brune    SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.
13899bd66eb0SPeter Brune 
13909bd66eb0SPeter Brune    Input Parameters:
13919bd66eb0SPeter Brune .  snes - nonlinear context obtained from SNESCreate()
13929bd66eb0SPeter Brune 
13939bd66eb0SPeter Brune    Output Parameters:
13949bd66eb0SPeter Brune +  projectfunc - function for projecting the function to the bounds
13959bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
13969bd66eb0SPeter Brune 
13979bd66eb0SPeter Brune    Logically Collective on SNES
13989bd66eb0SPeter Brune 
13999bd66eb0SPeter Brune     Level: developer
14009bd66eb0SPeter Brune 
14019bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, get, line search
14029bd66eb0SPeter Brune 
1403f1c6b773SPeter Brune .seealso: SNESLineSearchSetVIFunctions(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
14049bd66eb0SPeter Brune @*/
1405f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc)
14069bd66eb0SPeter Brune {
14079bd66eb0SPeter Brune   PetscFunctionBegin;
14089bd66eb0SPeter Brune   if (projectfunc) *projectfunc = linesearch->ops->viproject;
14099bd66eb0SPeter Brune   if (normfunc) *normfunc = linesearch->ops->vinorm;
14109bd66eb0SPeter Brune   PetscFunctionReturn(0);
14119bd66eb0SPeter Brune }
14129bd66eb0SPeter Brune 
14139bd66eb0SPeter Brune #undef __FUNCT__
1414f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchRegister"
1415bf7f4e0aSPeter Brune /*@C
1416f1c6b773SPeter Brune   SNESLineSearchRegister - See SNESLineSearchRegisterDynamic()
1417bf7f4e0aSPeter Brune 
1418bf7f4e0aSPeter Brune   Level: advanced
1419bf7f4e0aSPeter Brune @*/
1420f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNESLineSearch))
1421bf7f4e0aSPeter Brune {
1422bf7f4e0aSPeter Brune   char           fullname[PETSC_MAX_PATH_LEN];
1423bf7f4e0aSPeter Brune   PetscErrorCode ierr;
1424bf7f4e0aSPeter Brune 
1425bf7f4e0aSPeter Brune   PetscFunctionBegin;
1426bf7f4e0aSPeter Brune   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
1427f1c6b773SPeter Brune   ierr = PetscFListAdd(&SNESLineSearchList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
1428bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1429bf7f4e0aSPeter Brune }
1430