xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision 71eef1ae3c0a4c903a6a6abce6f6ad17aba1e92a)
1b45d2f2cSJed Brown #include <petsc-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 /*@
12cd7522eaSPeter Brune    SNESLineSearchCreate - Creates the line search context.
13f40b411bSPeter Brune 
14cd7522eaSPeter Brune    Logically Collective on Comm
15f40b411bSPeter Brune 
16f40b411bSPeter Brune    Input Parameters:
17cd7522eaSPeter Brune .  comm - MPI communicator for the line search (typically from the associated SNES context).
18f40b411bSPeter Brune 
19f40b411bSPeter Brune    Output Parameters:
208e557f58SPeter Brune .  outlinesearch - the new linesearch context
21f40b411bSPeter Brune 
2278bcb3b5SPeter Brune    Level: beginner
23f40b411bSPeter Brune 
24cd7522eaSPeter 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;
6259405d5eSPeter Brune   linesearch->max_its       = 1;
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 /*@
7178bcb3b5SPeter Brune    SNESLineSearchSetUp - Prepares the line search for being applied by allocating
7278bcb3b5SPeter Brune    any required vectors.
73f40b411bSPeter Brune 
74cd7522eaSPeter Brune    Collective on SNESLineSearch
75f40b411bSPeter Brune 
76f40b411bSPeter Brune    Input Parameters:
77f40b411bSPeter Brune .  linesearch - The LineSearch instance.
78f40b411bSPeter Brune 
79cd7522eaSPeter Brune    Notes:
80cd7522eaSPeter Brune    For most cases, this needn't be called outside of SNESLineSearchApply().
81cd7522eaSPeter Brune    The only current case where this is called outside of this is for the VI
8278bcb3b5SPeter Brune    solvers, which modify the solution and work vectors before the first call
83cd7522eaSPeter Brune    of SNESLineSearchApply, requiring the SNESLineSearch work vectors to be
84cd7522eaSPeter Brune    allocated upfront.
85cd7522eaSPeter Brune 
86cd7522eaSPeter Brune 
8778bcb3b5SPeter Brune    Level: advanced
88f40b411bSPeter Brune 
89f1c6b773SPeter Brune .keywords: SNESLineSearch, SetUp
90f40b411bSPeter Brune 
91f1c6b773SPeter Brune .seealso: SNESLineSearchReset()
92f40b411bSPeter Brune @*/
93f40b411bSPeter Brune 
94f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch) {
95bf7f4e0aSPeter Brune   PetscErrorCode ierr;
96bf7f4e0aSPeter Brune   PetscFunctionBegin;
97bf7f4e0aSPeter Brune   if (!((PetscObject)linesearch)->type_name) {
981a4f838cSPeter Brune     ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr);
99bf7f4e0aSPeter Brune   }
100bf7f4e0aSPeter Brune   if (!linesearch->setupcalled) {
1019bd66eb0SPeter Brune     if (!linesearch->vec_sol_new) {
102bf7f4e0aSPeter Brune       ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new);CHKERRQ(ierr);
1039bd66eb0SPeter Brune     }
1049bd66eb0SPeter Brune     if (!linesearch->vec_func_new) {
1059bd66eb0SPeter Brune       ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new);CHKERRQ(ierr);
1069bd66eb0SPeter Brune     }
107bf7f4e0aSPeter Brune     if (linesearch->ops->setup) {
108bf7f4e0aSPeter Brune       ierr = (*linesearch->ops->setup)(linesearch);CHKERRQ(ierr);
109bf7f4e0aSPeter Brune     }
110bf7f4e0aSPeter Brune     linesearch->lambda = linesearch->damping;
111bf7f4e0aSPeter Brune     linesearch->setupcalled = PETSC_TRUE;
112bf7f4e0aSPeter Brune   }
113bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
114bf7f4e0aSPeter Brune }
115bf7f4e0aSPeter Brune 
116bf7f4e0aSPeter Brune #undef __FUNCT__
117f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchReset"
118f40b411bSPeter Brune 
119f40b411bSPeter Brune /*@
12078bcb3b5SPeter Brune    SNESLineSearchReset - Undoes the SetUp and deletes any Vecs or Mats allocated by the line search.
121f40b411bSPeter Brune 
122f1c6b773SPeter Brune    Collective on SNESLineSearch
123f40b411bSPeter Brune 
124f40b411bSPeter Brune    Input Parameters:
125f40b411bSPeter Brune .  linesearch - The LineSearch instance.
126f40b411bSPeter Brune 
12778bcb3b5SPeter Brune    Level: intermediate
128f40b411bSPeter Brune 
129cd7522eaSPeter Brune .keywords: SNESLineSearch, Reset
130f40b411bSPeter Brune 
131f1c6b773SPeter Brune .seealso: SNESLineSearchSetUp()
132f40b411bSPeter Brune @*/
133f40b411bSPeter Brune 
134f1c6b773SPeter Brune PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch) {
135bf7f4e0aSPeter Brune   PetscErrorCode ierr;
136bf7f4e0aSPeter Brune   PetscFunctionBegin;
137bf7f4e0aSPeter Brune   if (linesearch->ops->reset) {
138bf7f4e0aSPeter Brune     (*linesearch->ops->reset)(linesearch);
139bf7f4e0aSPeter Brune   }
140bf7f4e0aSPeter Brune   ierr = VecDestroy(&linesearch->vec_sol_new);CHKERRQ(ierr);
141bf7f4e0aSPeter Brune   ierr = VecDestroy(&linesearch->vec_func_new);CHKERRQ(ierr);
142bf7f4e0aSPeter Brune 
143bf7f4e0aSPeter Brune   ierr = VecDestroyVecs(linesearch->nwork, &linesearch->work);CHKERRQ(ierr);
144bf7f4e0aSPeter Brune   linesearch->nwork = 0;
145bf7f4e0aSPeter Brune   linesearch->setupcalled = PETSC_FALSE;
146bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
147bf7f4e0aSPeter Brune }
148bf7f4e0aSPeter Brune 
14986d74e61SPeter Brune 
15086d74e61SPeter Brune #undef __FUNCT__
151f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetPreCheck"
15286d74e61SPeter Brune /*@C
153f1c6b773SPeter Brune    SNESLineSearchSetPreCheck - Sets a pre-check function for the line search routine.
15486d74e61SPeter Brune 
155f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
15686d74e61SPeter Brune 
15786d74e61SPeter Brune    Input Parameters:
158f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
15986d74e61SPeter Brune .  func       - [optional] function evaluation routine
16086d74e61SPeter Brune -  ctx        - [optional] user-defined context for private data for the
16186d74e61SPeter Brune                 function evaluation routine (may be PETSC_NULL)
16286d74e61SPeter Brune 
16386d74e61SPeter Brune    Calling sequence of func:
164f1c6b773SPeter Brune $    func (SNESLineSearch snes,Vec x,Vec y, PetscBool *changed);
16586d74e61SPeter Brune 
16686d74e61SPeter Brune +  x - solution vector
16786d74e61SPeter Brune .  y - search direction vector
168cd7522eaSPeter Brune -  changed - flag to indicate the precheck changed x or y.
16986d74e61SPeter Brune 
17086d74e61SPeter Brune    Level: intermediate
17186d74e61SPeter Brune 
17286d74e61SPeter Brune .keywords: set, linesearch, pre-check
17386d74e61SPeter Brune 
174f1c6b773SPeter Brune .seealso: SNESLineSearchSetPostCheck()
17586d74e61SPeter Brune @*/
176f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetPreCheck(SNESLineSearch linesearch, SNESLineSearchPreCheckFunc func,void *ctx)
17786d74e61SPeter Brune {
1789bd66eb0SPeter Brune   PetscFunctionBegin;
179f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
18086d74e61SPeter Brune   if (func) linesearch->ops->precheckstep = func;
18186d74e61SPeter Brune   if (ctx) linesearch->precheckctx = ctx;
18286d74e61SPeter Brune   PetscFunctionReturn(0);
18386d74e61SPeter Brune }
18486d74e61SPeter Brune 
18586d74e61SPeter Brune 
18686d74e61SPeter Brune #undef __FUNCT__
187f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetPreCheck"
18886d74e61SPeter Brune /*@C
189cd7522eaSPeter Brune    SNESLineSearchSetPreCheck - Gets the pre-check function for the line search routine.
19086d74e61SPeter Brune 
19186d74e61SPeter Brune    Input Parameters:
192f1c6b773SPeter Brune .  linesearch - the SNESLineSearch context
19386d74e61SPeter Brune 
19486d74e61SPeter Brune    Output Parameters:
19586d74e61SPeter Brune +  func       - [optional] function evaluation routine
19686d74e61SPeter Brune -  ctx        - [optional] user-defined context for private data for the
19786d74e61SPeter Brune                 function evaluation routine (may be PETSC_NULL)
19886d74e61SPeter Brune 
19986d74e61SPeter Brune    Level: intermediate
20086d74e61SPeter Brune 
20186d74e61SPeter Brune .keywords: get, linesearch, pre-check
20286d74e61SPeter Brune 
203f1c6b773SPeter Brune .seealso: SNESLineSearchGetPostCheck(), SNESLineSearchSetPreCheck()
20486d74e61SPeter Brune @*/
205f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetPreCheck(SNESLineSearch linesearch, SNESLineSearchPreCheckFunc *func,void **ctx)
20686d74e61SPeter Brune {
2079bd66eb0SPeter Brune   PetscFunctionBegin;
208f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
20986d74e61SPeter Brune   if (func) *func = linesearch->ops->precheckstep;
21086d74e61SPeter Brune   if (ctx) *ctx = linesearch->precheckctx;
21186d74e61SPeter Brune   PetscFunctionReturn(0);
21286d74e61SPeter Brune }
21386d74e61SPeter Brune 
21486d74e61SPeter Brune 
21586d74e61SPeter Brune #undef __FUNCT__
216f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetPostCheck"
21786d74e61SPeter Brune /*@C
218f1c6b773SPeter Brune    SNESLineSearchSetPostCheck - Sets a post-check function for the line search routine.
21986d74e61SPeter Brune 
220f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
22186d74e61SPeter Brune 
22286d74e61SPeter Brune    Input Parameters:
223f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
22486d74e61SPeter Brune .  func       - [optional] function evaluation routine
22586d74e61SPeter Brune -  ctx        - [optional] user-defined context for private data for the
22686d74e61SPeter Brune                 function evaluation routine (may be PETSC_NULL)
22786d74e61SPeter Brune 
22886d74e61SPeter Brune    Calling sequence of func:
229f1c6b773SPeter Brune $    func (SNESLineSearch linesearch,Vec x, Vec w, Vec y, PetscBool *changed_w, *changed_y);
23086d74e61SPeter Brune 
23186d74e61SPeter Brune +  x - old solution vector
23286d74e61SPeter Brune .  y - search direction vector
23386d74e61SPeter Brune .  w - new solution vector
2348e557f58SPeter Brune .  changed_y - indicates that the line search changed y
2358e557f58SPeter Brune .  changed_w - indicates that the line search changed w
23686d74e61SPeter Brune 
23786d74e61SPeter Brune    Level: intermediate
23886d74e61SPeter Brune 
23986d74e61SPeter Brune .keywords: set, linesearch, post-check
24086d74e61SPeter Brune 
241f1c6b773SPeter Brune .seealso: SNESLineSearchSetPreCheck()
24286d74e61SPeter Brune @*/
243f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetPostCheck(SNESLineSearch linesearch, SNESLineSearchPostCheckFunc func,void *ctx)
24486d74e61SPeter Brune {
24586d74e61SPeter Brune   PetscFunctionBegin;
246f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
24786d74e61SPeter Brune   if (func) linesearch->ops->postcheckstep = func;
24886d74e61SPeter Brune   if (ctx) linesearch->postcheckctx = ctx;
24986d74e61SPeter Brune   PetscFunctionReturn(0);
25086d74e61SPeter Brune }
25186d74e61SPeter Brune 
25286d74e61SPeter Brune 
25386d74e61SPeter Brune #undef __FUNCT__
254f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetPostCheck"
25586d74e61SPeter Brune /*@C
256f1c6b773SPeter Brune    SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine.
25786d74e61SPeter Brune 
25886d74e61SPeter Brune    Input Parameters:
259f1c6b773SPeter Brune .  linesearch - the SNESLineSearch context
26086d74e61SPeter Brune 
26186d74e61SPeter Brune    Output Parameters:
26286d74e61SPeter Brune +  func       - [optional] function evaluation routine
26386d74e61SPeter Brune -  ctx        - [optional] user-defined context for private data for the
26486d74e61SPeter Brune                 function evaluation routine (may be PETSC_NULL)
26586d74e61SPeter Brune 
26686d74e61SPeter Brune    Level: intermediate
26786d74e61SPeter Brune 
26886d74e61SPeter Brune .keywords: get, linesearch, post-check
26986d74e61SPeter Brune 
270f1c6b773SPeter Brune .seealso: SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck()
27186d74e61SPeter Brune @*/
272f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetPostCheck(SNESLineSearch linesearch, SNESLineSearchPostCheckFunc *func,void **ctx)
27386d74e61SPeter Brune {
2749bd66eb0SPeter Brune   PetscFunctionBegin;
275f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
27686d74e61SPeter Brune   if (func) *func = linesearch->ops->postcheckstep;
27786d74e61SPeter Brune   if (ctx) *ctx = linesearch->postcheckctx;
27886d74e61SPeter Brune   PetscFunctionReturn(0);
27986d74e61SPeter Brune }
28086d74e61SPeter Brune 
28186d74e61SPeter Brune 
282bf7f4e0aSPeter Brune #undef __FUNCT__
283f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPreCheck"
284f40b411bSPeter Brune /*@
285f1c6b773SPeter Brune    SNESLineSearchPreCheck - Prepares the line search for being applied.
286f40b411bSPeter Brune 
287cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
288f40b411bSPeter Brune 
289f40b411bSPeter Brune    Input Parameters:
290f40b411bSPeter Brune .  linesearch - The linesearch instance.
291f40b411bSPeter Brune 
292f40b411bSPeter Brune    Output Parameters:
2938e557f58SPeter Brune .  changed - Indicator that the precheck routine has changed anything
294f40b411bSPeter Brune 
295f40b411bSPeter Brune    Level: Beginner
296f40b411bSPeter Brune 
297f1c6b773SPeter Brune .keywords: SNESLineSearch, Create
298f40b411bSPeter Brune 
299f1c6b773SPeter Brune .seealso: SNESLineSearchPostCheck()
300f40b411bSPeter Brune @*/
301f1c6b773SPeter Brune PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch, PetscBool * changed)
302bf7f4e0aSPeter Brune {
303bf7f4e0aSPeter Brune   PetscErrorCode ierr;
304bf7f4e0aSPeter Brune   PetscFunctionBegin;
305bf7f4e0aSPeter Brune   *changed = PETSC_FALSE;
306bf7f4e0aSPeter Brune   if (linesearch->ops->precheckstep) {
307c87759e9SPeter Brune     ierr = (*linesearch->ops->precheckstep)(linesearch, linesearch->vec_sol, linesearch->vec_update, changed, linesearch->precheckctx);CHKERRQ(ierr);
308bf7f4e0aSPeter Brune   }
309bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
310bf7f4e0aSPeter Brune }
311bf7f4e0aSPeter Brune 
312bf7f4e0aSPeter Brune #undef __FUNCT__
313f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPostCheck"
314f40b411bSPeter Brune /*@
315f1c6b773SPeter Brune    SNESLineSearchPostCheck - Prepares the line search for being applied.
316f40b411bSPeter Brune 
317cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
318f40b411bSPeter Brune 
319f40b411bSPeter Brune    Input Parameters:
3208e557f58SPeter Brune .  linesearch - The linesearch context
321f40b411bSPeter Brune 
322f40b411bSPeter Brune    Output Parameters:
32378bcb3b5SPeter Brune +  changed_Y - Indicator if the direction Y has been changed.
32478bcb3b5SPeter Brune -  changed_W - Indicator if the new candidate solution W has been changed.
325f40b411bSPeter Brune 
326f40b411bSPeter Brune    Level: Intermediate
327f40b411bSPeter Brune 
328f1c6b773SPeter Brune .keywords: SNESLineSearch, Create
329f40b411bSPeter Brune 
330f1c6b773SPeter Brune .seealso: SNESLineSearchPreCheck()
331f40b411bSPeter Brune @*/
332f1c6b773SPeter Brune PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch, PetscBool * changed_Y, PetscBool * changed_W)
333bf7f4e0aSPeter Brune {
334bf7f4e0aSPeter Brune   PetscErrorCode ierr;
335bf7f4e0aSPeter Brune   PetscFunctionBegin;
336bf7f4e0aSPeter Brune   *changed_Y = PETSC_FALSE;
337bf7f4e0aSPeter Brune   *changed_W = PETSC_FALSE;
338bf7f4e0aSPeter Brune   if (linesearch->ops->postcheckstep) {
339c87759e9SPeter Brune     ierr = (*linesearch->ops->postcheckstep)(linesearch, linesearch->vec_sol, linesearch->vec_update, linesearch->vec_sol_new, changed_Y, changed_W, linesearch->postcheckctx);CHKERRQ(ierr);
34086d74e61SPeter Brune   }
34186d74e61SPeter Brune   PetscFunctionReturn(0);
34286d74e61SPeter Brune }
34386d74e61SPeter Brune 
34486d74e61SPeter Brune 
34586d74e61SPeter Brune #undef __FUNCT__
346f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPreCheckPicard"
34786d74e61SPeter Brune /*@C
34886d74e61SPeter Brune    SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration
34986d74e61SPeter Brune 
350cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
35186d74e61SPeter Brune 
35286d74e61SPeter Brune    Input Arguments:
35386d74e61SPeter Brune +  linesearch - linesearch context
35486d74e61SPeter Brune .  X - base state for this step
35586d74e61SPeter Brune .  Y - initial correction
35686d74e61SPeter Brune 
35786d74e61SPeter Brune    Output Arguments:
35886d74e61SPeter Brune +  Y - correction, possibly modified
35986d74e61SPeter Brune -  changed - flag indicating that Y was modified
36086d74e61SPeter Brune 
36186d74e61SPeter Brune    Options Database Key:
362cd7522eaSPeter Brune +  -snes_linesearch_precheck_picard - activate this routine
363cd7522eaSPeter Brune -  -snes_linesearch_precheck_picard_angle - angle
36486d74e61SPeter Brune 
36586d74e61SPeter Brune    Level: advanced
36686d74e61SPeter Brune 
36786d74e61SPeter Brune    Notes:
36886d74e61SPeter Brune    This function should be passed to SNESLineSearchSetPreCheck()
36986d74e61SPeter Brune 
37086d74e61SPeter Brune    The justification for this method involves the linear convergence of a Picard iteration
37186d74e61SPeter Brune    so the Picard linearization should be provided in place of the "Jacobian". This correction
37286d74e61SPeter Brune    is generally not useful when using a Newton linearization.
37386d74e61SPeter Brune 
37486d74e61SPeter Brune    Reference:
37586d74e61SPeter Brune    Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology.
37686d74e61SPeter Brune 
37786d74e61SPeter Brune .seealso: SNESLineSearchSetPreCheck()
37886d74e61SPeter Brune @*/
379f1c6b773SPeter Brune PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed,void *ctx)
38086d74e61SPeter Brune {
38186d74e61SPeter Brune   PetscErrorCode ierr;
38286d74e61SPeter Brune   PetscReal      angle = *(PetscReal*)linesearch->precheckctx;
38386d74e61SPeter Brune   Vec            Ylast;
38486d74e61SPeter Brune   PetscScalar    dot;
385c87759e9SPeter Brune 
38686d74e61SPeter Brune   PetscInt       iter;
38786d74e61SPeter Brune   PetscReal      ynorm,ylastnorm,theta,angle_radians;
38886d74e61SPeter Brune   SNES           snes;
38986d74e61SPeter Brune 
39086d74e61SPeter Brune   PetscFunctionBegin;
391f1c6b773SPeter Brune   ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
39286d74e61SPeter Brune   ierr = PetscObjectQuery((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject*)&Ylast);CHKERRQ(ierr);
39386d74e61SPeter Brune   if (!Ylast) {
39486d74e61SPeter Brune     ierr = VecDuplicate(Y,&Ylast);CHKERRQ(ierr);
39586d74e61SPeter Brune     ierr = PetscObjectCompose((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject)Ylast);CHKERRQ(ierr);
39686d74e61SPeter Brune     ierr = PetscObjectDereference((PetscObject)Ylast);CHKERRQ(ierr);
39786d74e61SPeter Brune   }
39886d74e61SPeter Brune   ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr);
39986d74e61SPeter Brune   if (iter < 2) {
40086d74e61SPeter Brune     ierr = VecCopy(Y,Ylast);CHKERRQ(ierr);
40186d74e61SPeter Brune     *changed = PETSC_FALSE;
40286d74e61SPeter Brune     PetscFunctionReturn(0);
40386d74e61SPeter Brune   }
40486d74e61SPeter Brune 
40586d74e61SPeter Brune   ierr = VecDot(Y,Ylast,&dot);CHKERRQ(ierr);
40686d74e61SPeter Brune   ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr);
40786d74e61SPeter Brune   ierr = VecNorm(Ylast,NORM_2,&ylastnorm);CHKERRQ(ierr);
40886d74e61SPeter Brune   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
40986d74e61SPeter Brune   theta = acos((double)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0));
41086d74e61SPeter Brune   angle_radians = angle * PETSC_PI / 180.;
41186d74e61SPeter Brune   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
41286d74e61SPeter Brune     /* Modify the step Y */
41386d74e61SPeter Brune     PetscReal alpha,ydiffnorm;
41486d74e61SPeter Brune     ierr = VecAXPY(Ylast,-1.0,Y);CHKERRQ(ierr);
41586d74e61SPeter Brune     ierr = VecNorm(Ylast,NORM_2,&ydiffnorm);CHKERRQ(ierr);
41686d74e61SPeter Brune     alpha = ylastnorm / ydiffnorm;
41786d74e61SPeter Brune     ierr = VecCopy(Y,Ylast);CHKERRQ(ierr);
41886d74e61SPeter Brune     ierr = VecScale(Y,alpha);CHKERRQ(ierr);
41986d74e61SPeter 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);
42086d74e61SPeter Brune   } else {
42186d74e61SPeter Brune     ierr = PetscInfo2(snes,"Angle %G degrees exceeds threshold %G, no correction applied\n",theta*180./PETSC_PI,angle);CHKERRQ(ierr);
42286d74e61SPeter Brune     ierr = VecCopy(Y,Ylast);CHKERRQ(ierr);
42386d74e61SPeter Brune     *changed = PETSC_FALSE;
424bf7f4e0aSPeter Brune   }
425bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
426bf7f4e0aSPeter Brune }
427bf7f4e0aSPeter Brune 
428bf7f4e0aSPeter Brune #undef __FUNCT__
429f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchApply"
430f40b411bSPeter Brune /*@
431cd7522eaSPeter Brune    SNESLineSearchApply - Computes the line-search update.
432f40b411bSPeter Brune 
433f1c6b773SPeter Brune    Collective on SNESLineSearch
434f40b411bSPeter Brune 
435f40b411bSPeter Brune    Input Parameters:
4368e557f58SPeter Brune +  linesearch - The linesearch context
4378e557f58SPeter Brune .  X - The current solution
4388e557f58SPeter Brune .  F - The current function
4398e557f58SPeter Brune .  fnorm - The current norm
4408e557f58SPeter Brune -  Y - The search direction
441f40b411bSPeter Brune 
442f40b411bSPeter Brune    Output Parameters:
4438e557f58SPeter Brune +  X - The new solution
4448e557f58SPeter Brune .  F - The new function
4458e557f58SPeter Brune -  fnorm - The new function norm
446f40b411bSPeter Brune 
447cd7522eaSPeter Brune    Options Database Keys:
448cd7522eaSPeter Brune + -snes_linesearch_type - The Line search method
449cd7522eaSPeter Brune . -snes_linesearch_monitor - Print progress of line searches
450cd7522eaSPeter Brune . -snes_linesearch_damping - The linesearch damping parameter.
451cd7522eaSPeter Brune . -snes_linesearch_norms   - Turn on/off the linesearch norms
452cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess.
453cd7522eaSPeter Brune - -snes_linesearch_max_it - The number of iterations for iterative line searches.
454cd7522eaSPeter Brune 
455cd7522eaSPeter Brune    Notes:
456cd7522eaSPeter Brune    This is typically called from within a SNESSolve() implementation in order to
457cd7522eaSPeter Brune    help with convergence of the nonlinear method.  Various SNES types use line searches
458cd7522eaSPeter Brune    in different ways, but the overarching theme is that a line search is used to determine
459cd7522eaSPeter Brune    an optimal damping parameter of a step at each iteration of the method.  Each
460cd7522eaSPeter Brune    application of the line search may invoke SNESComputeFunction several times, and
461cd7522eaSPeter Brune    therefore may be fairly expensive.
462cd7522eaSPeter Brune 
463f40b411bSPeter Brune    Level: Intermediate
464f40b411bSPeter Brune 
465f1c6b773SPeter Brune .keywords: SNESLineSearch, Create
466f40b411bSPeter Brune 
467cd7522eaSPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchPreCheck(), SNESLineSearchPostCheck(), SNESSolve(), SNESComputeFunction()
468f40b411bSPeter Brune @*/
469f1c6b773SPeter Brune PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y) {
470bf7f4e0aSPeter Brune   PetscErrorCode ierr;
471bf7f4e0aSPeter Brune   PetscFunctionBegin;
472bf7f4e0aSPeter Brune 
473bf7f4e0aSPeter Brune   /* check the pointers */
474f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
475bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
476bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
477bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
478bf7f4e0aSPeter Brune 
479bf7f4e0aSPeter Brune   linesearch->success = PETSC_TRUE;
480bf7f4e0aSPeter Brune 
481bf7f4e0aSPeter Brune   linesearch->vec_sol = X;
482bf7f4e0aSPeter Brune   linesearch->vec_update = Y;
483bf7f4e0aSPeter Brune   linesearch->vec_func = F;
484bf7f4e0aSPeter Brune 
485f1c6b773SPeter Brune   ierr = SNESLineSearchSetUp(linesearch);CHKERRQ(ierr);
486bf7f4e0aSPeter Brune 
487bf7f4e0aSPeter Brune   if (!linesearch->keeplambda)
488bf7f4e0aSPeter Brune     linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */
489bf7f4e0aSPeter Brune 
490bf7f4e0aSPeter Brune   if (fnorm) {
491bf7f4e0aSPeter Brune     linesearch->fnorm = *fnorm;
492bf7f4e0aSPeter Brune   } else {
493bf7f4e0aSPeter Brune     ierr = VecNorm(F, NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
494bf7f4e0aSPeter Brune   }
495bf7f4e0aSPeter Brune 
496f1c6b773SPeter Brune   ierr = PetscLogEventBegin(SNESLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
497bf7f4e0aSPeter Brune 
498bf7f4e0aSPeter Brune   ierr = (*linesearch->ops->apply)(linesearch);CHKERRQ(ierr);
499bf7f4e0aSPeter Brune 
500f1c6b773SPeter Brune   ierr = PetscLogEventEnd(SNESLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
501bf7f4e0aSPeter Brune 
502bf7f4e0aSPeter Brune   if (fnorm)
503bf7f4e0aSPeter Brune     *fnorm = linesearch->fnorm;
504bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
505bf7f4e0aSPeter Brune }
506bf7f4e0aSPeter Brune 
507bf7f4e0aSPeter Brune #undef __FUNCT__
508f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchDestroy"
509f40b411bSPeter Brune /*@
510f1c6b773SPeter Brune    SNESLineSearchDestroy - Destroys the line search instance.
511f40b411bSPeter Brune 
512f1c6b773SPeter Brune    Collective on SNESLineSearch
513f40b411bSPeter Brune 
514f40b411bSPeter Brune    Input Parameters:
5158e557f58SPeter Brune .  linesearch - The linesearch context
516f40b411bSPeter Brune 
517f40b411bSPeter Brune    Level: Intermediate
518f40b411bSPeter Brune 
51978bcb3b5SPeter Brune .keywords: SNESLineSearch, Destroy
520f40b411bSPeter Brune 
521cd7522eaSPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchReset(), SNESDestroy()
522f40b411bSPeter Brune @*/
523f1c6b773SPeter Brune PetscErrorCode SNESLineSearchDestroy(SNESLineSearch * linesearch) {
524bf7f4e0aSPeter Brune   PetscErrorCode ierr;
525bf7f4e0aSPeter Brune   PetscFunctionBegin;
526bf7f4e0aSPeter Brune   if (!*linesearch) PetscFunctionReturn(0);
527f1c6b773SPeter Brune   PetscValidHeaderSpecific((*linesearch),SNESLINESEARCH_CLASSID,1);
528bf7f4e0aSPeter Brune   if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = 0; PetscFunctionReturn(0);}
529bf7f4e0aSPeter Brune   ierr = PetscObjectDepublish((*linesearch));CHKERRQ(ierr);
530f1c6b773SPeter Brune   ierr = SNESLineSearchReset(*linesearch);
531bf7f4e0aSPeter Brune   if ((*linesearch)->ops->destroy) {
532bf7f4e0aSPeter Brune     (*linesearch)->ops->destroy(*linesearch);
533bf7f4e0aSPeter Brune   }
534bf7f4e0aSPeter Brune   ierr = PetscViewerDestroy(&(*linesearch)->monitor);CHKERRQ(ierr);
535e7058c64SPeter Brune   ierr = PetscHeaderDestroy(linesearch);CHKERRQ(ierr);
536bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
537bf7f4e0aSPeter Brune }
538bf7f4e0aSPeter Brune 
539bf7f4e0aSPeter Brune #undef __FUNCT__
540f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetMonitor"
541f40b411bSPeter Brune /*@
542cd7522eaSPeter Brune    SNESLineSearchSetMonitor - Turns on/off printing useful information and debugging output about the line search.
543bf7f4e0aSPeter Brune 
544bf7f4e0aSPeter Brune    Input Parameters:
545bf7f4e0aSPeter Brune +  snes - nonlinear context obtained from SNESCreate()
546bf7f4e0aSPeter Brune -  flg - PETSC_TRUE to monitor the line search
547bf7f4e0aSPeter Brune 
548bf7f4e0aSPeter Brune    Logically Collective on SNES
549bf7f4e0aSPeter Brune 
550bf7f4e0aSPeter Brune    Options Database:
5518e557f58SPeter Brune .   -snes_linesearch_monitor - enables the monitor
552bf7f4e0aSPeter Brune 
553bf7f4e0aSPeter Brune    Level: intermediate
554bf7f4e0aSPeter Brune 
555bf7f4e0aSPeter Brune 
556cd7522eaSPeter Brune .seealso: SNESLineSearchGetMonitor(), PetscViewer
557bf7f4e0aSPeter Brune @*/
558f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetMonitor(SNESLineSearch linesearch, PetscBool flg)
559bf7f4e0aSPeter Brune {
560bf7f4e0aSPeter Brune 
561bf7f4e0aSPeter Brune   PetscErrorCode ierr;
562bf7f4e0aSPeter Brune   PetscFunctionBegin;
563bf7f4e0aSPeter Brune   if (flg && !linesearch->monitor) {
564bf7f4e0aSPeter Brune     ierr = PetscViewerASCIIOpen(((PetscObject)linesearch)->comm,"stdout",&linesearch->monitor);CHKERRQ(ierr);
565bf7f4e0aSPeter Brune   } else if (!flg && linesearch->monitor) {
566bf7f4e0aSPeter Brune     ierr = PetscViewerDestroy(&linesearch->monitor);CHKERRQ(ierr);
567bf7f4e0aSPeter Brune   }
568bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
569bf7f4e0aSPeter Brune }
570bf7f4e0aSPeter Brune 
571bf7f4e0aSPeter Brune #undef __FUNCT__
572f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetMonitor"
573f40b411bSPeter Brune /*@
574cd7522eaSPeter Brune    SNESLineSearchGetMonitor - Gets the PetscViewer instance for the line search monitor.
5756a388c36SPeter Brune 
576f40b411bSPeter Brune    Input Parameters:
5778e557f58SPeter Brune .  linesearch - linesearch context
578f40b411bSPeter Brune 
579f40b411bSPeter Brune    Input Parameters:
5808e557f58SPeter Brune .  monitor - monitor context
581f40b411bSPeter Brune 
582f40b411bSPeter Brune    Logically Collective on SNES
583f40b411bSPeter Brune 
584f40b411bSPeter Brune 
585f40b411bSPeter Brune    Options Database Keys:
5868e557f58SPeter Brune .   -snes_linesearch_monitor - enables the monitor
587f40b411bSPeter Brune 
588f40b411bSPeter Brune    Level: intermediate
589f40b411bSPeter Brune 
590f40b411bSPeter Brune 
591cd7522eaSPeter Brune .seealso: SNESLineSearchSetMonitor(), PetscViewer
592f40b411bSPeter Brune @*/
593f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetMonitor(SNESLineSearch linesearch, PetscViewer *monitor)
5946a388c36SPeter Brune {
5956a388c36SPeter Brune 
5966a388c36SPeter Brune   PetscFunctionBegin;
597f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
5986a388c36SPeter Brune   if (monitor) {
5996a388c36SPeter Brune     PetscValidPointer(monitor, 2);
6006a388c36SPeter Brune     *monitor = linesearch->monitor;
6016a388c36SPeter Brune   }
6026a388c36SPeter Brune   PetscFunctionReturn(0);
6036a388c36SPeter Brune }
6046a388c36SPeter Brune 
6056a388c36SPeter Brune #undef __FUNCT__
606f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetFromOptions"
607f40b411bSPeter Brune /*@
608f1c6b773SPeter Brune    SNESLineSearchSetFromOptions - Sets options for the line search
609f40b411bSPeter Brune 
610f40b411bSPeter Brune    Input Parameters:
6118e557f58SPeter Brune .  linesearch - linesearch context
612f40b411bSPeter Brune 
613f40b411bSPeter Brune    Options Database Keys:
614cd7522eaSPeter Brune + -snes_linesearch_type - The Line search method
61571eef1aeSPeter Brune . -snes_linesearch_minlambda - The minimum step length
616*1a9b3a06SPeter Brune . -snes_linesearch_maxstep - The maximum step size
617*1a9b3a06SPeter Brune . -snes_linesearch_rtol - Relative tolerance for iterative line searches
618*1a9b3a06SPeter Brune . -snes_linesearch_atol - Absolute tolerance for iterative line searches
619*1a9b3a06SPeter Brune . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches
620*1a9b3a06SPeter Brune . -snes_linesearch_max_it - The number of iterations for iterative line searches
621cd7522eaSPeter Brune . -snes_linesearch_monitor - Print progress of line searches
6228e557f58SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter
623cd7522eaSPeter Brune . -snes_linesearch_norms   - Turn on/off the linesearch norms
624cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess.
62559405d5eSPeter Brune . -snes_linesearch_order - The polynomial order of approximation used in the linesearch
626*1a9b3a06SPeter Brune . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method
627*1a9b3a06SPeter Brune - -snes_linesearch_precheck_picard_angle - Angle used in picard precheck method
628f40b411bSPeter Brune 
629f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
630f40b411bSPeter Brune 
631f40b411bSPeter Brune    Level: intermediate
632f40b411bSPeter Brune 
633f40b411bSPeter Brune 
634f1c6b773SPeter Brune .seealso: SNESLineSearchCreate()
635f40b411bSPeter Brune @*/
636f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch) {
637bf7f4e0aSPeter Brune   PetscErrorCode ierr;
63859405d5eSPeter Brune   const char       *orders[] = {"linear", "quadratic", "cubic"};
63959405d5eSPeter Brune   PetscInt         indx = 0;
6401a4f838cSPeter Brune   const char     *deft = SNESLINESEARCHBASIC;
641bf7f4e0aSPeter Brune   char           type[256];
642bf7f4e0aSPeter Brune   PetscBool      flg, set;
643bf7f4e0aSPeter Brune   PetscFunctionBegin;
644f1c6b773SPeter Brune   if (!SNESLineSearchRegisterAllCalled) {ierr = SNESLineSearchRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
645bf7f4e0aSPeter Brune 
646bf7f4e0aSPeter Brune   ierr = PetscObjectOptionsBegin((PetscObject)linesearch);CHKERRQ(ierr);
647bf7f4e0aSPeter Brune   if (((PetscObject)linesearch)->type_name) {
648bf7f4e0aSPeter Brune     deft = ((PetscObject)linesearch)->type_name;
649bf7f4e0aSPeter Brune   }
6507a35526eSPeter Brune   ierr = PetscOptionsList("-snes_linesearch_type","Line-search method","SNESLineSearchSetType",SNESLineSearchList,deft,type,256,&flg);CHKERRQ(ierr);
651bf7f4e0aSPeter Brune   if (flg) {
652f1c6b773SPeter Brune     ierr = SNESLineSearchSetType(linesearch,type);CHKERRQ(ierr);
653bf7f4e0aSPeter Brune   } else if (!((PetscObject)linesearch)->type_name) {
654f1c6b773SPeter Brune     ierr = SNESLineSearchSetType(linesearch,deft);CHKERRQ(ierr);
655bf7f4e0aSPeter Brune   }
656bf7f4e0aSPeter Brune   if (linesearch->ops->setfromoptions) {
657bf7f4e0aSPeter Brune     (*linesearch->ops->setfromoptions)(linesearch);CHKERRQ(ierr);
658bf7f4e0aSPeter Brune   }
659bf7f4e0aSPeter Brune 
6607a35526eSPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_monitor","Print progress of line searches","SNESSNESLineSearchSetMonitor",
661bf7f4e0aSPeter Brune                           linesearch->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
662f1c6b773SPeter Brune   if (set) {ierr = SNESLineSearchSetMonitor(linesearch,flg);CHKERRQ(ierr);}
663bf7f4e0aSPeter Brune 
664*1a9b3a06SPeter Brune   /* tolerances */
66571eef1aeSPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_minlambda","Minimum step length","SNESLineSearchSetTolerances",linesearch->steptol,&linesearch->steptol,0);CHKERRQ(ierr);
666*1a9b3a06SPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_maxstep","Maximum step size","SNESLineSearchSetTolerances",linesearch->maxstep,&linesearch->maxstep,0);CHKERRQ(ierr);
667*1a9b3a06SPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_rtol","Relative tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->rtol,&linesearch->rtol,0);CHKERRQ(ierr);
668*1a9b3a06SPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_atol","Absolute tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->atol,&linesearch->atol,0);CHKERRQ(ierr);
669*1a9b3a06SPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_ltol","Change in lambda tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->ltol,&linesearch->ltol,0);CHKERRQ(ierr);
6708e557f58SPeter Brune   ierr = PetscOptionsInt("-snes_linesearch_max_it","Maximum iterations for iterative line searches","SNESLineSearchSetTolerances",linesearch->max_its,&linesearch->max_its,0);CHKERRQ(ierr);
6717a35526eSPeter Brune 
672*1a9b3a06SPeter Brune   /* damping parameters */
673*1a9b3a06SPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_damping","Line search damping and initial step guess","SNESLineSearchSetDamping",linesearch->damping,&linesearch->damping,0);CHKERRQ(ierr);
674*1a9b3a06SPeter Brune 
675*1a9b3a06SPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_keeplambda","Use previous lambda as damping","SNESLineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,0);CHKERRQ(ierr);
676*1a9b3a06SPeter Brune 
677*1a9b3a06SPeter Brune   /* precheck */
6787a35526eSPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_precheck_picard","Use a correction that sometimes improves convergence of Picard iteration","SNESLineSearchPreCheckPicard",flg,&flg,&set);CHKERRQ(ierr);
6797a35526eSPeter Brune   if (set) {
6807a35526eSPeter Brune     if (flg) {
6817a35526eSPeter Brune       linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */
6827a35526eSPeter Brune       ierr = PetscOptionsReal("-snes_linesearch_precheck_picard_angle","Maximum angle at which to activate the correction",
6837a35526eSPeter Brune                               "none",linesearch->precheck_picard_angle,&linesearch->precheck_picard_angle,PETSC_NULL);CHKERRQ(ierr);
6847a35526eSPeter Brune       ierr = SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&linesearch->precheck_picard_angle);CHKERRQ(ierr);
6857a35526eSPeter Brune     } else {
6867a35526eSPeter Brune       ierr = SNESLineSearchSetPreCheck(linesearch,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
6877a35526eSPeter Brune     }
6887a35526eSPeter Brune   }
689*1a9b3a06SPeter Brune 
6908e557f58SPeter Brune   ierr = PetscOptionsEList("-snes_linesearch_order",  "Order of approximation", "SNESLineSearchSetOrder", orders,3,orders[(PetscInt)linesearch->order],&indx,&flg);CHKERRQ(ierr);
69159405d5eSPeter Brune   if (flg) {
69259405d5eSPeter Brune     switch (indx) {
69359405d5eSPeter Brune     case 0: linesearch->order = SNES_LINESEARCH_LINEAR;
69459405d5eSPeter Brune       break;
69559405d5eSPeter Brune     case 1: linesearch->order = SNES_LINESEARCH_QUADRATIC;
69659405d5eSPeter Brune       break;
69759405d5eSPeter Brune     case 2: linesearch->order = SNES_LINESEARCH_CUBIC;
69859405d5eSPeter Brune       break;
69959405d5eSPeter Brune     }
70059405d5eSPeter Brune   }
70159405d5eSPeter Brune 
702*1a9b3a06SPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_norms","Compute final norms in line search","SNESLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,0);CHKERRQ(ierr);
7037a35526eSPeter Brune 
704bf7f4e0aSPeter Brune   ierr = PetscObjectProcessOptionsHandlers((PetscObject)linesearch);CHKERRQ(ierr);
705bf7f4e0aSPeter Brune   ierr = PetscOptionsEnd();CHKERRQ(ierr);
706bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
707bf7f4e0aSPeter Brune }
708bf7f4e0aSPeter Brune 
709bf7f4e0aSPeter Brune #undef __FUNCT__
710f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchView"
711f40b411bSPeter Brune /*@
712cd7522eaSPeter Brune    SNESLineSearchView - Prints useful information about the line search not
713cd7522eaSPeter Brune    related to an individual call.
714f40b411bSPeter Brune 
715f40b411bSPeter Brune    Input Parameters:
7168e557f58SPeter Brune .  linesearch - linesearch context
717f40b411bSPeter Brune 
718f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
719f40b411bSPeter Brune 
720f40b411bSPeter Brune    Level: intermediate
721f40b411bSPeter Brune 
722f1c6b773SPeter Brune .seealso: SNESLineSearchCreate()
723f40b411bSPeter Brune @*/
7247f1410a3SPeter Brune PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer) {
7257f1410a3SPeter Brune   PetscErrorCode ierr;
7267f1410a3SPeter Brune   PetscBool      iascii;
727bf7f4e0aSPeter Brune   PetscFunctionBegin;
7287f1410a3SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
7297f1410a3SPeter Brune   if (!viewer) {
7307f1410a3SPeter Brune     ierr = PetscViewerASCIIGetStdout(((PetscObject)linesearch)->comm,&viewer);CHKERRQ(ierr);
7317f1410a3SPeter Brune   }
7327f1410a3SPeter Brune   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
7337f1410a3SPeter Brune   PetscCheckSameComm(linesearch,1,viewer,2);
734f40b411bSPeter Brune 
7357f1410a3SPeter Brune   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
7367f1410a3SPeter Brune   if (iascii) {
7377f1410a3SPeter Brune     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)linesearch,viewer,"SNESLineSearch Object");CHKERRQ(ierr);
7387f1410a3SPeter Brune     if (linesearch->ops->view) {
7397f1410a3SPeter Brune       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
7407f1410a3SPeter Brune       ierr = (*linesearch->ops->view)(linesearch,viewer);CHKERRQ(ierr);
7417f1410a3SPeter Brune       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
7427f1410a3SPeter Brune     }
7437f1410a3SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  maxstep=%G, minlambda=%G\n", linesearch->maxstep,linesearch->steptol);CHKERRQ(ierr);
7447f1410a3SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%G, absolute=%G, lambda=%G\n", linesearch->rtol,linesearch->atol,linesearch->ltol);CHKERRQ(ierr);
7457f1410a3SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D\n", linesearch->max_its);CHKERRQ(ierr);
7467f1410a3SPeter Brune     if (linesearch->ops->precheckstep) {
7477f1410a3SPeter Brune       if (linesearch->ops->precheckstep == SNESLineSearchPreCheckPicard) {
7487f1410a3SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  using precheck step to speed up Picard convergence\n", linesearch->max_its);CHKERRQ(ierr);
7497f1410a3SPeter Brune       } else {
7507f1410a3SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  using user-defined precheck step\n", linesearch->max_its);CHKERRQ(ierr);
7517f1410a3SPeter Brune       }
7527f1410a3SPeter Brune     }
7537f1410a3SPeter Brune     if (linesearch->ops->postcheckstep) {
7547f1410a3SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer,"  using user-defined postcheck step\n", linesearch->max_its);CHKERRQ(ierr);
7557f1410a3SPeter Brune     }
7567f1410a3SPeter Brune   }
757bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
758bf7f4e0aSPeter Brune }
759bf7f4e0aSPeter Brune 
760bf7f4e0aSPeter Brune #undef __FUNCT__
761f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetType"
762ea5d4fccSPeter Brune /*@C
763f1c6b773SPeter Brune    SNESLineSearchSetType - Sets the linesearch type
764f40b411bSPeter Brune 
765f40b411bSPeter Brune    Input Parameters:
7668e557f58SPeter Brune +  linesearch - linesearch context
767f40b411bSPeter Brune -  type - The type of line search to be used
768f40b411bSPeter Brune 
769cd7522eaSPeter Brune    Available Types:
770cd7522eaSPeter Brune +  basic - Simple damping line search.
7718e557f58SPeter Brune .  bt - Backtracking line search over the L2 norm of the function
7728e557f58SPeter Brune .  l2 - Secant line search over the L2 norm of the function
7738e557f58SPeter Brune .  cp - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x)
7748e557f58SPeter Brune -  shell - User provided SNESLineSearch implementation
775cd7522eaSPeter Brune 
776f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
777f40b411bSPeter Brune 
778f40b411bSPeter Brune    Level: intermediate
779f40b411bSPeter Brune 
780f40b411bSPeter Brune 
781f1c6b773SPeter Brune .seealso: SNESLineSearchCreate()
782f40b411bSPeter Brune @*/
783f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, const SNESLineSearchType type)
784bf7f4e0aSPeter Brune {
785bf7f4e0aSPeter Brune 
786f1c6b773SPeter Brune   PetscErrorCode ierr,(*r)(SNESLineSearch);
787bf7f4e0aSPeter Brune   PetscBool      match;
788bf7f4e0aSPeter Brune 
789bf7f4e0aSPeter Brune   PetscFunctionBegin;
790f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
791bf7f4e0aSPeter Brune   PetscValidCharPointer(type,2);
792bf7f4e0aSPeter Brune 
793bf7f4e0aSPeter Brune   ierr = PetscTypeCompare((PetscObject)linesearch,type,&match);CHKERRQ(ierr);
794bf7f4e0aSPeter Brune   if (match) PetscFunctionReturn(0);
795bf7f4e0aSPeter Brune 
796f1c6b773SPeter Brune   ierr =  PetscFListFind(SNESLineSearchList,((PetscObject)linesearch)->comm,type,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr);
797bf7f4e0aSPeter Brune   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type);
798bf7f4e0aSPeter Brune   /* Destroy the previous private linesearch context */
799bf7f4e0aSPeter Brune   if (linesearch->ops->destroy) {
800bf7f4e0aSPeter Brune     ierr = (*(linesearch)->ops->destroy)(linesearch);CHKERRQ(ierr);
801bf7f4e0aSPeter Brune     linesearch->ops->destroy = PETSC_NULL;
802bf7f4e0aSPeter Brune   }
803f1c6b773SPeter Brune   /* Reinitialize function pointers in SNESLineSearchOps structure */
804bf7f4e0aSPeter Brune   linesearch->ops->apply          = 0;
805bf7f4e0aSPeter Brune   linesearch->ops->view           = 0;
806bf7f4e0aSPeter Brune   linesearch->ops->setfromoptions = 0;
807bf7f4e0aSPeter Brune   linesearch->ops->destroy        = 0;
808bf7f4e0aSPeter Brune 
809bf7f4e0aSPeter Brune   ierr = PetscObjectChangeTypeName((PetscObject)linesearch,type);CHKERRQ(ierr);
810bf7f4e0aSPeter Brune   ierr = (*r)(linesearch);CHKERRQ(ierr);
811bf7f4e0aSPeter Brune #if defined(PETSC_HAVE_AMS)
812bf7f4e0aSPeter Brune   if (PetscAMSPublishAll) {
813bf7f4e0aSPeter Brune     ierr = PetscObjectAMSPublish((PetscObject)linesearch);CHKERRQ(ierr);
814bf7f4e0aSPeter Brune   }
815bf7f4e0aSPeter Brune #endif
816bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
817bf7f4e0aSPeter Brune }
818bf7f4e0aSPeter Brune 
819bf7f4e0aSPeter Brune #undef __FUNCT__
820f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetSNES"
821f40b411bSPeter Brune /*@
82278bcb3b5SPeter Brune    SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation.
823f40b411bSPeter Brune 
824f40b411bSPeter Brune    Input Parameters:
8258e557f58SPeter Brune +  linesearch - linesearch context
826f40b411bSPeter Brune -  snes - The snes instance
827f40b411bSPeter Brune 
82878bcb3b5SPeter Brune    Level: developer
82978bcb3b5SPeter Brune 
83078bcb3b5SPeter Brune    Notes:
83178bcb3b5SPeter Brune    This happens automatically when the line search is gotten/created with
83278bcb3b5SPeter Brune    SNESGetSNESLineSearch().  This routine is therefore mainly called within SNES
83378bcb3b5SPeter Brune    implementations.
834f40b411bSPeter Brune 
835f40b411bSPeter Brune 
836cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
837f40b411bSPeter Brune @*/
838f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes){
839bf7f4e0aSPeter Brune   PetscFunctionBegin;
840f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
841bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,2);
842bf7f4e0aSPeter Brune   linesearch->snes = snes;
843bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
844bf7f4e0aSPeter Brune }
845bf7f4e0aSPeter Brune 
846bf7f4e0aSPeter Brune #undef __FUNCT__
847f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetSNES"
848f40b411bSPeter Brune /*@
84978bcb3b5SPeter Brune    SNESLineSearchGetSNES - Gets the SNES for the linesearch for function evaluation.
850f40b411bSPeter Brune 
851f40b411bSPeter Brune    Input Parameters:
8528e557f58SPeter Brune .  linesearch - linesearch context
853f40b411bSPeter Brune 
854f40b411bSPeter Brune    Output Parameters:
855f40b411bSPeter Brune .  snes - The snes instance
856f40b411bSPeter Brune 
85778bcb3b5SPeter Brune    Level: advanced
858f40b411bSPeter Brune 
859cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
860f40b411bSPeter Brune @*/
861f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes){
862bf7f4e0aSPeter Brune   PetscFunctionBegin;
863f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
8646a388c36SPeter Brune   PetscValidPointer(snes, 2);
865bf7f4e0aSPeter Brune   *snes = linesearch->snes;
866bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
867bf7f4e0aSPeter Brune }
868bf7f4e0aSPeter Brune 
8696a388c36SPeter Brune #undef __FUNCT__
870f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetLambda"
871f40b411bSPeter Brune /*@
872f1c6b773SPeter Brune    SNESLineSearchGetLambda - Gets the last linesearch steplength discovered.
873f40b411bSPeter Brune 
874f40b411bSPeter Brune    Input Parameters:
8758e557f58SPeter Brune .  linesearch - linesearch context
876f40b411bSPeter Brune 
877f40b411bSPeter Brune    Output Parameters:
878cd7522eaSPeter Brune .  lambda - The last steplength computed during SNESLineSearchApply()
879f40b411bSPeter Brune 
88078bcb3b5SPeter Brune    Level: advanced
88178bcb3b5SPeter Brune 
8828e557f58SPeter Brune    Notes:
8838e557f58SPeter Brune    This is useful in methods where the solver is ill-scaled and
88478bcb3b5SPeter Brune    requires some adaptive notion of the difference in scale between the
88578bcb3b5SPeter Brune    solution and the function.  For instance, SNESQN may be scaled by the
88678bcb3b5SPeter Brune    line search lambda using the argument -snes_qn_scaling ls.
88778bcb3b5SPeter Brune 
888f40b411bSPeter Brune 
889cd7522eaSPeter Brune .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetDamping(), SNESLineSearchApply()
890f40b411bSPeter Brune @*/
891f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda)
8926a388c36SPeter Brune {
8936a388c36SPeter Brune   PetscFunctionBegin;
894f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
8956a388c36SPeter Brune   PetscValidPointer(lambda, 2);
8966a388c36SPeter Brune   *lambda = linesearch->lambda;
8976a388c36SPeter Brune   PetscFunctionReturn(0);
8986a388c36SPeter Brune }
8996a388c36SPeter Brune 
9006a388c36SPeter Brune #undef __FUNCT__
901f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetLambda"
902f40b411bSPeter Brune /*@
903f1c6b773SPeter Brune    SNESLineSearchSetLambda - Sets the linesearch steplength.
904f40b411bSPeter Brune 
905f40b411bSPeter Brune    Input Parameters:
9068e557f58SPeter Brune +  linesearch - linesearch context
907f40b411bSPeter Brune -  lambda - The last steplength.
908f40b411bSPeter Brune 
909cd7522eaSPeter Brune    Notes:
910cd7522eaSPeter Brune    This routine is typically used within implementations of SNESLineSearchApply
911cd7522eaSPeter Brune    to set the final steplength.  This routine (and SNESLineSearchGetLambda()) were
912cd7522eaSPeter Brune    added in order to facilitate Quasi-Newton methods that use the previous steplength
913cd7522eaSPeter Brune    as an inner scaling parameter.
914cd7522eaSPeter Brune 
91578bcb3b5SPeter Brune    Level: advanced
916f40b411bSPeter Brune 
917f1c6b773SPeter Brune .seealso: SNESLineSearchGetLambda()
918f40b411bSPeter Brune @*/
919f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
9206a388c36SPeter Brune {
9216a388c36SPeter Brune   PetscFunctionBegin;
922f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
9236a388c36SPeter Brune   linesearch->lambda = lambda;
9246a388c36SPeter Brune   PetscFunctionReturn(0);
9256a388c36SPeter Brune }
9266a388c36SPeter Brune 
9276a388c36SPeter Brune #undef  __FUNCT__
928f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetTolerances"
929f40b411bSPeter Brune /*@
93078bcb3b5SPeter Brune    SNESLineSearchGetTolerances - Gets the tolerances for the method.  These include
93178bcb3b5SPeter Brune    tolerances for the relative and absolute change in the function norm, the change
93278bcb3b5SPeter Brune    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
93378bcb3b5SPeter Brune    and the maximum number of iterations the line search procedure may take.
934f40b411bSPeter Brune 
935f40b411bSPeter Brune    Input Parameters:
9368e557f58SPeter Brune .  linesearch - linesearch context
937f40b411bSPeter Brune 
938f40b411bSPeter Brune    Output Parameters:
939516fe3c3SPeter Brune +  steptol - The minimum steplength
9406cc8e53bSPeter Brune .  maxstep - The maximum steplength
941516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
942516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
943516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
944516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
945f40b411bSPeter Brune 
94678bcb3b5SPeter Brune    Level: intermediate
94778bcb3b5SPeter Brune 
94878bcb3b5SPeter Brune    Notes:
94978bcb3b5SPeter Brune    Different line searches may implement these parameters slightly differently as
95078bcb3b5SPeter Brune    the method requires.
951516fe3c3SPeter Brune 
952f1c6b773SPeter Brune .seealso: SNESLineSearchSetTolerances()
953f40b411bSPeter Brune @*/
954f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
9556a388c36SPeter Brune {
9566a388c36SPeter Brune   PetscFunctionBegin;
957f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
958516fe3c3SPeter Brune   if (steptol) {
9596a388c36SPeter Brune     PetscValidPointer(steptol, 2);
9606a388c36SPeter Brune     *steptol = linesearch->steptol;
961516fe3c3SPeter Brune   }
962516fe3c3SPeter Brune   if (maxstep) {
963516fe3c3SPeter Brune     PetscValidPointer(maxstep, 3);
964516fe3c3SPeter Brune     *maxstep = linesearch->maxstep;
965516fe3c3SPeter Brune   }
966516fe3c3SPeter Brune   if (rtol) {
967516fe3c3SPeter Brune     PetscValidPointer(rtol, 4);
968516fe3c3SPeter Brune     *rtol = linesearch->rtol;
969516fe3c3SPeter Brune   }
970516fe3c3SPeter Brune   if (atol) {
971516fe3c3SPeter Brune     PetscValidPointer(atol, 5);
972516fe3c3SPeter Brune     *atol = linesearch->atol;
973516fe3c3SPeter Brune   }
974516fe3c3SPeter Brune   if (ltol) {
975516fe3c3SPeter Brune     PetscValidPointer(ltol, 6);
976516fe3c3SPeter Brune     *ltol = linesearch->ltol;
977516fe3c3SPeter Brune   }
978516fe3c3SPeter Brune   if (max_its) {
979516fe3c3SPeter Brune     PetscValidPointer(max_its, 7);
980516fe3c3SPeter Brune     *max_its = linesearch->max_its;
981516fe3c3SPeter Brune   }
9826a388c36SPeter Brune   PetscFunctionReturn(0);
9836a388c36SPeter Brune }
9846a388c36SPeter Brune 
9856a388c36SPeter Brune #undef  __FUNCT__
986f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetTolerances"
987f40b411bSPeter Brune /*@
98878bcb3b5SPeter Brune    SNESLineSearchSetTolerances -  Gets the tolerances for the method.  These include
98978bcb3b5SPeter Brune    tolerances for the relative and absolute change in the function norm, the change
99078bcb3b5SPeter Brune    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
99178bcb3b5SPeter Brune    and the maximum number of iterations the line search procedure may take.
992f40b411bSPeter Brune 
993f40b411bSPeter Brune    Input Parameters:
9948e557f58SPeter Brune +  linesearch - linesearch context
995516fe3c3SPeter Brune .  steptol - The minimum steplength
9966cc8e53bSPeter Brune .  maxstep - The maximum steplength
997516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
998516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
999516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
1000516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
1001f40b411bSPeter Brune 
100278bcb3b5SPeter Brune    Notes:
100378bcb3b5SPeter Brune    The user may choose to not set any of the tolerances in the method using PETSC_DEFAULT in
100478bcb3b5SPeter Brune    place of an argument.
1005f40b411bSPeter Brune 
100678bcb3b5SPeter Brune    Level: intermediate
1007516fe3c3SPeter Brune 
1008f1c6b773SPeter Brune .seealso: SNESLineSearchGetTolerances()
1009f40b411bSPeter Brune @*/
1010f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its)
10116a388c36SPeter Brune {
10126a388c36SPeter Brune   PetscFunctionBegin;
1013f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1014d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,steptol,2);
1015d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,maxstep,3);
1016d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,rtol,4);
1017d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,atol,5);
1018d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,ltol,6);
1019d3952184SSatish Balay   PetscValidLogicalCollectiveInt(linesearch,max_its,7);
1020d3952184SSatish Balay 
1021d3952184SSatish Balay   if ( steptol!= PETSC_DEFAULT) {
1022d3952184SSatish Balay     if ( steptol < 0.0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Minimum step length %G must be non-negative",steptol);
10236a388c36SPeter Brune     linesearch->steptol = steptol;
1024d3952184SSatish Balay   }
1025d3952184SSatish Balay 
1026d3952184SSatish Balay   if ( maxstep!= PETSC_DEFAULT) {
1027d3952184SSatish Balay     if ( maxstep < 0.0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Maximum step length %G must be non-negative",maxstep);
1028516fe3c3SPeter Brune     linesearch->maxstep = maxstep;
1029d3952184SSatish Balay   }
1030d3952184SSatish Balay 
1031d3952184SSatish Balay   if (rtol != PETSC_DEFAULT) {
1032d3952184SSatish Balay     if (rtol < 0.0 || 1.0 <= rtol) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Relative tolerance %G must be non-negative and less than 1.0",rtol);
1033516fe3c3SPeter Brune     linesearch->rtol = rtol;
1034d3952184SSatish Balay   }
1035d3952184SSatish Balay 
1036d3952184SSatish Balay   if (atol != PETSC_DEFAULT) {
1037d3952184SSatish Balay     if (atol < 0.0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %G must be non-negative",atol);
1038516fe3c3SPeter Brune     linesearch->atol = atol;
1039d3952184SSatish Balay   }
1040d3952184SSatish Balay 
1041d3952184SSatish Balay   if (ltol != PETSC_DEFAULT) {
1042d3952184SSatish Balay     if (ltol < 0.0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Labmda tolerance %G must be non-negative",ltol);
1043516fe3c3SPeter Brune   linesearch->ltol = ltol;
1044d3952184SSatish Balay   }
1045d3952184SSatish Balay 
1046d3952184SSatish Balay   if (max_its != PETSC_DEFAULT) {
1047d3952184SSatish Balay     if (max_its < 0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",max_its);
1048516fe3c3SPeter Brune     linesearch->max_its = max_its;
1049d3952184SSatish Balay   }
1050d3952184SSatish Balay 
10516a388c36SPeter Brune   PetscFunctionReturn(0);
10526a388c36SPeter Brune }
10536a388c36SPeter Brune 
1054516fe3c3SPeter Brune 
10556a388c36SPeter Brune #undef __FUNCT__
1056f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetDamping"
1057f40b411bSPeter Brune /*@
1058f1c6b773SPeter Brune    SNESLineSearchGetDamping - Gets the line search damping parameter.
1059f40b411bSPeter Brune 
1060f40b411bSPeter Brune    Input Parameters:
10618e557f58SPeter Brune .  linesearch - linesearch context
1062f40b411bSPeter Brune 
1063f40b411bSPeter Brune    Output Parameters:
10648e557f58SPeter Brune .  damping - The damping parameter
1065f40b411bSPeter Brune 
106678bcb3b5SPeter Brune    Level: advanced
1067f40b411bSPeter Brune 
106878bcb3b5SPeter Brune .seealso: SNESLineSearchGetStepTolerance(), SNESQN
1069f40b411bSPeter Brune @*/
1070f40b411bSPeter Brune 
1071f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping)
10726a388c36SPeter Brune {
10736a388c36SPeter Brune   PetscFunctionBegin;
1074f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
10756a388c36SPeter Brune   PetscValidPointer(damping, 2);
10766a388c36SPeter Brune   *damping = linesearch->damping;
10776a388c36SPeter Brune   PetscFunctionReturn(0);
10786a388c36SPeter Brune }
10796a388c36SPeter Brune 
10806a388c36SPeter Brune #undef __FUNCT__
1081f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetDamping"
1082f40b411bSPeter Brune /*@
1083f1c6b773SPeter Brune    SNESLineSearchSetDamping - Sets the line search damping paramter.
1084f40b411bSPeter Brune 
1085f40b411bSPeter Brune    Input Parameters:
108678bcb3b5SPeter Brune .  linesearch - linesearch context
108778bcb3b5SPeter Brune .  damping - The damping parameter
1088f40b411bSPeter Brune 
1089f40b411bSPeter Brune    Level: intermediate
1090f40b411bSPeter Brune 
1091cd7522eaSPeter Brune    Notes:
1092cd7522eaSPeter Brune    The basic line search merely takes the update step scaled by the damping parameter.
1093cd7522eaSPeter Brune    The use of the damping parameter in the l2 and cp line searches is much more subtle;
109478bcb3b5SPeter Brune    it is used as a starting point in calculating the secant step. However, the eventual
1095cd7522eaSPeter Brune    step may be of greater length than the damping parameter.  In the bt line search it is
1096cd7522eaSPeter Brune    used as the maximum possible step length, as the bt line search only backtracks.
1097cd7522eaSPeter Brune 
1098f1c6b773SPeter Brune .seealso: SNESLineSearchGetDamping()
1099f40b411bSPeter Brune @*/
1100f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping)
11016a388c36SPeter Brune {
11026a388c36SPeter Brune   PetscFunctionBegin;
1103f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
11046a388c36SPeter Brune   linesearch->damping = damping;
11056a388c36SPeter Brune   PetscFunctionReturn(0);
11066a388c36SPeter Brune }
11076a388c36SPeter Brune 
11086a388c36SPeter Brune #undef __FUNCT__
110959405d5eSPeter Brune #define __FUNCT__ "SNESLineSearchGetOrder"
111059405d5eSPeter Brune /*@
111159405d5eSPeter Brune    SNESLineSearchGetOrder - Gets the line search approximation order.
111259405d5eSPeter Brune 
111359405d5eSPeter Brune    Input Parameters:
111478bcb3b5SPeter Brune .  linesearch - linesearch context
111559405d5eSPeter Brune 
111659405d5eSPeter Brune    Output Parameters:
11178e557f58SPeter Brune .  order - The order
111859405d5eSPeter Brune 
111978bcb3b5SPeter Brune    Possible Values for order:
11208e557f58SPeter Brune +  SNES_LINESEARCH_LINEAR - linear method
11218e557f58SPeter Brune .  SNES_LINESEARCH_QUADRATIC - quadratic method
11228e557f58SPeter Brune -  SNES_LINESEARCH_CUBIC - cubic method
112378bcb3b5SPeter Brune 
112459405d5eSPeter Brune    Level: intermediate
112559405d5eSPeter Brune 
112659405d5eSPeter Brune .seealso: SNESLineSearchSetOrder()
112759405d5eSPeter Brune @*/
112859405d5eSPeter Brune 
112959405d5eSPeter Brune PetscErrorCode  SNESLineSearchGetOrder(SNESLineSearch linesearch,SNESLineSearchOrder *order)
113059405d5eSPeter Brune {
113159405d5eSPeter Brune   PetscFunctionBegin;
113259405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
113359405d5eSPeter Brune   PetscValidPointer(order, 2);
113459405d5eSPeter Brune   *order = linesearch->order;
113559405d5eSPeter Brune   PetscFunctionReturn(0);
113659405d5eSPeter Brune }
113759405d5eSPeter Brune 
113859405d5eSPeter Brune #undef __FUNCT__
113959405d5eSPeter Brune #define __FUNCT__ "SNESLineSearchSetOrder"
114059405d5eSPeter Brune /*@
114159405d5eSPeter Brune    SNESLineSearchSetOrder - Sets the line search damping paramter.
114259405d5eSPeter Brune 
114359405d5eSPeter Brune    Input Parameters:
114478bcb3b5SPeter Brune .  linesearch - linesearch context
114578bcb3b5SPeter Brune .  order - The damping parameter
114659405d5eSPeter Brune 
114759405d5eSPeter Brune    Level: intermediate
114859405d5eSPeter Brune 
114978bcb3b5SPeter Brune    Possible Values for order:
11508e557f58SPeter Brune +  SNES_LINESEARCH_LINEAR - linear method
11518e557f58SPeter Brune .  SNES_LINESEARCH_QUADRATIC - quadratic method
11528e557f58SPeter Brune -  SNES_LINESEARCH_CUBIC - cubic method
115378bcb3b5SPeter Brune 
115459405d5eSPeter Brune    Notes:
115559405d5eSPeter Brune    Variable orders are supported by the following line searches:
115678bcb3b5SPeter Brune +  bt - cubic and quadratic
115778bcb3b5SPeter Brune -  cp - linear and quadratic
115859405d5eSPeter Brune 
115959405d5eSPeter Brune .seealso: SNESLineSearchGetOrder()
116059405d5eSPeter Brune @*/
116159405d5eSPeter Brune PetscErrorCode  SNESLineSearchSetOrder(SNESLineSearch linesearch,SNESLineSearchOrder order)
116259405d5eSPeter Brune {
116359405d5eSPeter Brune   PetscFunctionBegin;
116459405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
116559405d5eSPeter Brune   linesearch->order = order;
116659405d5eSPeter Brune   PetscFunctionReturn(0);
116759405d5eSPeter Brune }
116859405d5eSPeter Brune 
116959405d5eSPeter Brune #undef __FUNCT__
1170f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetNorms"
1171f40b411bSPeter Brune /*@
1172f1c6b773SPeter Brune    SNESLineSearchGetNorms - Gets the norms for for X, Y, and F.
1173f40b411bSPeter Brune 
1174f40b411bSPeter Brune    Input Parameters:
117578bcb3b5SPeter Brune .  linesearch - linesearch context
1176f40b411bSPeter Brune 
1177f40b411bSPeter Brune    Output Parameters:
1178f40b411bSPeter Brune +  xnorm - The norm of the current solution
1179f40b411bSPeter Brune .  fnorm - The norm of the current function
1180f40b411bSPeter Brune -  ynorm - The norm of the current update
1181f40b411bSPeter Brune 
1182cd7522eaSPeter Brune    Notes:
1183cd7522eaSPeter Brune    This function is mainly called from SNES implementations.
1184cd7522eaSPeter Brune 
118578bcb3b5SPeter Brune    Level: developer
1186f40b411bSPeter Brune 
1187f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms() SNESLineSearchGetVecs()
1188f40b411bSPeter Brune @*/
1189f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm)
1190bf7f4e0aSPeter Brune {
1191bf7f4e0aSPeter Brune   PetscFunctionBegin;
1192f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1193bf7f4e0aSPeter Brune   if (xnorm) {
1194bf7f4e0aSPeter Brune     *xnorm = linesearch->xnorm;
1195bf7f4e0aSPeter Brune   }
1196bf7f4e0aSPeter Brune   if (fnorm) {
1197bf7f4e0aSPeter Brune     *fnorm = linesearch->fnorm;
1198bf7f4e0aSPeter Brune   }
1199bf7f4e0aSPeter Brune   if (ynorm) {
1200bf7f4e0aSPeter Brune     *ynorm = linesearch->ynorm;
1201bf7f4e0aSPeter Brune   }
1202bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1203bf7f4e0aSPeter Brune }
1204bf7f4e0aSPeter Brune 
1205e7058c64SPeter Brune #undef __FUNCT__
1206f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetNorms"
1207f40b411bSPeter Brune /*@
1208f1c6b773SPeter Brune    SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F.
1209f40b411bSPeter Brune 
1210f40b411bSPeter Brune    Input Parameters:
121178bcb3b5SPeter Brune +  linesearch - linesearch context
1212f40b411bSPeter Brune .  xnorm - The norm of the current solution
1213f40b411bSPeter Brune .  fnorm - The norm of the current function
1214f40b411bSPeter Brune -  ynorm - The norm of the current update
1215f40b411bSPeter Brune 
121678bcb3b5SPeter Brune    Level: advanced
1217f40b411bSPeter Brune 
1218f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1219f40b411bSPeter Brune @*/
1220f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
12216a388c36SPeter Brune {
12226a388c36SPeter Brune   PetscFunctionBegin;
1223f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
12246a388c36SPeter Brune   linesearch->xnorm = xnorm;
12256a388c36SPeter Brune   linesearch->fnorm = fnorm;
12266a388c36SPeter Brune   linesearch->ynorm = ynorm;
12276a388c36SPeter Brune   PetscFunctionReturn(0);
12286a388c36SPeter Brune }
12296a388c36SPeter Brune 
12306a388c36SPeter Brune #undef __FUNCT__
1231f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchComputeNorms"
1232f40b411bSPeter Brune /*@
1233f1c6b773SPeter Brune    SNESLineSearchComputeNorms - Computes the norms of X, F, and Y.
1234f40b411bSPeter Brune 
1235f40b411bSPeter Brune    Input Parameters:
123678bcb3b5SPeter Brune .  linesearch - linesearch context
1237f40b411bSPeter Brune 
1238f40b411bSPeter Brune    Options Database Keys:
12398e557f58SPeter Brune .   -snes_linesearch_norms - turn norm computation on or off
1240f40b411bSPeter Brune 
1241f40b411bSPeter Brune    Level: intermediate
1242f40b411bSPeter Brune 
124378bcb3b5SPeter Brune .seealso: SNESLineSearchGetNorms, SNESLineSearchSetNorms(), SNESLineSearchSetComputeNorms()
1244f40b411bSPeter Brune @*/
1245f1c6b773SPeter Brune PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
12466a388c36SPeter Brune {
12476a388c36SPeter Brune   PetscErrorCode ierr;
12489bd66eb0SPeter Brune   SNES snes;
12496a388c36SPeter Brune   PetscFunctionBegin;
12506a388c36SPeter Brune   if (linesearch->norms) {
12519bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
1252f1c6b773SPeter Brune       ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
12539bd66eb0SPeter Brune       ierr = VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
12549bd66eb0SPeter Brune       ierr = VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
12559bd66eb0SPeter Brune       ierr = (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);CHKERRQ(ierr);
12569bd66eb0SPeter Brune     } else {
12576a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_func,   NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
12586a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_sol,    NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
12596a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
12606a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_func,     NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
12616a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_sol,      NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
12626a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_update,   NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
12636a388c36SPeter Brune     }
12649bd66eb0SPeter Brune   }
12656a388c36SPeter Brune   PetscFunctionReturn(0);
12666a388c36SPeter Brune }
12676a388c36SPeter Brune 
12686f263ca3SPeter Brune 
12696f263ca3SPeter Brune #undef __FUNCT__
12706f263ca3SPeter Brune #define __FUNCT__ "SNESLineSearchSetComputeNorms"
12716f263ca3SPeter Brune /*@
12726f263ca3SPeter Brune    SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search.
12736f263ca3SPeter Brune 
12746f263ca3SPeter Brune    Input Parameters:
127578bcb3b5SPeter Brune +  linesearch  - linesearch context
127678bcb3b5SPeter Brune -  flg  - indicates whether or not to compute norms
12776f263ca3SPeter Brune 
12786f263ca3SPeter Brune    Options Database Keys:
12798e557f58SPeter Brune .   -snes_linesearch_norms - turn norm computation on or off
12806f263ca3SPeter Brune 
12816f263ca3SPeter Brune    Notes:
12821a4f838cSPeter Brune    This is most relevant to the SNESLINESEARCHBASIC line search type.
12836f263ca3SPeter Brune 
12846f263ca3SPeter Brune    Level: intermediate
12856f263ca3SPeter Brune 
12861a4f838cSPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetNorms(), SNESLineSearchComputeNorms(), SNESLINESEARCHBASIC
12876f263ca3SPeter Brune @*/
12886f263ca3SPeter Brune PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg)
12896f263ca3SPeter Brune {
12906f263ca3SPeter Brune   PetscFunctionBegin;
12916f263ca3SPeter Brune   linesearch->norms = flg;
12926f263ca3SPeter Brune   PetscFunctionReturn(0);
12936f263ca3SPeter Brune }
12946f263ca3SPeter Brune 
12956a388c36SPeter Brune #undef __FUNCT__
1296f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetVecs"
1297f40b411bSPeter Brune /*@
1298f1c6b773SPeter Brune    SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context
1299f40b411bSPeter Brune 
1300f40b411bSPeter Brune    Input Parameters:
130178bcb3b5SPeter Brune .  linesearch - linesearch context
1302f40b411bSPeter Brune 
1303f40b411bSPeter Brune    Output Parameters:
1304f40b411bSPeter Brune +  X - The old solution
1305f40b411bSPeter Brune .  F - The old function
1306f40b411bSPeter Brune .  Y - The search direction
1307f40b411bSPeter Brune .  W - The new solution
1308f40b411bSPeter Brune -  G - The new function
1309f40b411bSPeter Brune 
131078bcb3b5SPeter Brune    Level: advanced
1311f40b411bSPeter Brune 
1312f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1313f40b411bSPeter Brune @*/
1314f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G) {
13156a388c36SPeter Brune   PetscFunctionBegin;
1316f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
13176a388c36SPeter Brune   if (X) {
13186a388c36SPeter Brune     PetscValidPointer(X, 2);
13196a388c36SPeter Brune     *X = linesearch->vec_sol;
13206a388c36SPeter Brune   }
13216a388c36SPeter Brune   if (F) {
13226a388c36SPeter Brune     PetscValidPointer(F, 3);
13236a388c36SPeter Brune     *F = linesearch->vec_func;
13246a388c36SPeter Brune   }
13256a388c36SPeter Brune   if (Y) {
13266a388c36SPeter Brune     PetscValidPointer(Y, 4);
13276a388c36SPeter Brune     *Y = linesearch->vec_update;
13286a388c36SPeter Brune   }
13296a388c36SPeter Brune   if (W) {
13306a388c36SPeter Brune     PetscValidPointer(W, 5);
13316a388c36SPeter Brune     *W = linesearch->vec_sol_new;
13326a388c36SPeter Brune   }
13336a388c36SPeter Brune   if (G) {
13346a388c36SPeter Brune     PetscValidPointer(G, 6);
13356a388c36SPeter Brune     *G = linesearch->vec_func_new;
13366a388c36SPeter Brune   }
13376a388c36SPeter Brune 
13386a388c36SPeter Brune   PetscFunctionReturn(0);
13396a388c36SPeter Brune }
13406a388c36SPeter Brune 
13416a388c36SPeter Brune #undef __FUNCT__
1342f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetVecs"
1343f40b411bSPeter Brune /*@
1344f1c6b773SPeter Brune    SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context
1345f40b411bSPeter Brune 
1346f40b411bSPeter Brune    Input Parameters:
134778bcb3b5SPeter Brune +  linesearch - linesearch context
1348f40b411bSPeter Brune .  X - The old solution
1349f40b411bSPeter Brune .  F - The old function
1350f40b411bSPeter Brune .  Y - The search direction
1351f40b411bSPeter Brune .  W - The new solution
1352f40b411bSPeter Brune -  G - The new function
1353f40b411bSPeter Brune 
135478bcb3b5SPeter Brune    Level: advanced
1355f40b411bSPeter Brune 
1356f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms(), SNESLineSearchGetVecs()
1357f40b411bSPeter Brune @*/
1358f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G) {
13596a388c36SPeter Brune   PetscFunctionBegin;
1360f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
13616a388c36SPeter Brune   if (X) {
13626a388c36SPeter Brune     PetscValidHeaderSpecific(X,VEC_CLASSID,2);
13636a388c36SPeter Brune     linesearch->vec_sol = X;
13646a388c36SPeter Brune   }
13656a388c36SPeter Brune   if (F) {
13666a388c36SPeter Brune     PetscValidHeaderSpecific(F,VEC_CLASSID,3);
13676a388c36SPeter Brune     linesearch->vec_func = F;
13686a388c36SPeter Brune   }
13696a388c36SPeter Brune   if (Y) {
13706a388c36SPeter Brune     PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
13716a388c36SPeter Brune     linesearch->vec_update = Y;
13726a388c36SPeter Brune   }
13736a388c36SPeter Brune   if (W) {
13746a388c36SPeter Brune     PetscValidHeaderSpecific(W,VEC_CLASSID,5);
13756a388c36SPeter Brune     linesearch->vec_sol_new = W;
13766a388c36SPeter Brune   }
13776a388c36SPeter Brune   if (G) {
13786a388c36SPeter Brune     PetscValidHeaderSpecific(G,VEC_CLASSID,6);
13796a388c36SPeter Brune     linesearch->vec_func_new = G;
13806a388c36SPeter Brune   }
13816a388c36SPeter Brune 
13826a388c36SPeter Brune   PetscFunctionReturn(0);
13836a388c36SPeter Brune }
13846a388c36SPeter Brune 
13856a388c36SPeter Brune #undef __FUNCT__
1386f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchAppendOptionsPrefix"
1387e7058c64SPeter Brune /*@C
1388f1c6b773SPeter Brune    SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1389e7058c64SPeter Brune    SNES options in the database.
1390e7058c64SPeter Brune 
1391cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
1392e7058c64SPeter Brune 
1393e7058c64SPeter Brune    Input Parameters:
1394e7058c64SPeter Brune +  snes - the SNES context
1395e7058c64SPeter Brune -  prefix - the prefix to prepend to all option names
1396e7058c64SPeter Brune 
1397e7058c64SPeter Brune    Notes:
1398e7058c64SPeter Brune    A hyphen (-) must NOT be given at the beginning of the prefix name.
1399e7058c64SPeter Brune    The first character of all runtime options is AUTOMATICALLY the hyphen.
1400e7058c64SPeter Brune 
1401e7058c64SPeter Brune    Level: advanced
1402e7058c64SPeter Brune 
1403f1c6b773SPeter Brune .keywords: SNESLineSearch, append, options, prefix, database
1404e7058c64SPeter Brune 
1405e7058c64SPeter Brune .seealso: SNESGetOptionsPrefix()
1406e7058c64SPeter Brune @*/
1407f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[])
1408e7058c64SPeter Brune {
1409e7058c64SPeter Brune   PetscErrorCode ierr;
1410e7058c64SPeter Brune 
1411e7058c64SPeter Brune   PetscFunctionBegin;
1412f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1413e7058c64SPeter Brune   ierr = PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1414e7058c64SPeter Brune   PetscFunctionReturn(0);
1415e7058c64SPeter Brune }
1416e7058c64SPeter Brune 
1417e7058c64SPeter Brune #undef __FUNCT__
1418f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetOptionsPrefix"
1419e7058c64SPeter Brune /*@C
1420f1c6b773SPeter Brune    SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all
1421f1c6b773SPeter Brune    SNESLineSearch options in the database.
1422e7058c64SPeter Brune 
1423e7058c64SPeter Brune    Not Collective
1424e7058c64SPeter Brune 
1425e7058c64SPeter Brune    Input Parameter:
1426cd7522eaSPeter Brune .  linesearch - the SNESLineSearch context
1427e7058c64SPeter Brune 
1428e7058c64SPeter Brune    Output Parameter:
1429e7058c64SPeter Brune .  prefix - pointer to the prefix string used
1430e7058c64SPeter Brune 
14318e557f58SPeter Brune    Notes:
14328e557f58SPeter Brune    On the fortran side, the user should pass in a string 'prefix' of
1433e7058c64SPeter Brune    sufficient length to hold the prefix.
1434e7058c64SPeter Brune 
1435e7058c64SPeter Brune    Level: advanced
1436e7058c64SPeter Brune 
1437f1c6b773SPeter Brune .keywords: SNESLineSearch, get, options, prefix, database
1438e7058c64SPeter Brune 
1439e7058c64SPeter Brune .seealso: SNESAppendOptionsPrefix()
1440e7058c64SPeter Brune @*/
1441f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[])
1442e7058c64SPeter Brune {
1443e7058c64SPeter Brune   PetscErrorCode ierr;
1444e7058c64SPeter Brune 
1445e7058c64SPeter Brune   PetscFunctionBegin;
1446f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1447e7058c64SPeter Brune   ierr = PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1448e7058c64SPeter Brune   PetscFunctionReturn(0);
1449e7058c64SPeter Brune }
1450bf7f4e0aSPeter Brune 
1451bf7f4e0aSPeter Brune #undef __FUNCT__
1452f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetWork"
1453f40b411bSPeter Brune /*@
1454f1c6b773SPeter Brune    SNESLineSearchGetWork - Gets work vectors for the line search.
1455f40b411bSPeter Brune 
1456f40b411bSPeter Brune    Input Parameter:
1457f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
1458f40b411bSPeter Brune -  nwork - the number of work vectors
1459f40b411bSPeter Brune 
1460f40b411bSPeter Brune    Level: developer
1461f40b411bSPeter Brune 
1462cd7522eaSPeter Brune    Notes:
1463cd7522eaSPeter Brune    This is typically called at the beginning of a SNESLineSearch or SNESLineSearchShell implementation.
1464cd7522eaSPeter Brune 
1465f1c6b773SPeter Brune .keywords: SNESLineSearch, work, vector
1466f40b411bSPeter Brune 
1467f40b411bSPeter Brune .seealso: SNESDefaultGetWork()
1468f40b411bSPeter Brune @*/
1469f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetWork(SNESLineSearch linesearch, PetscInt nwork)
1470bf7f4e0aSPeter Brune {
1471bf7f4e0aSPeter Brune   PetscErrorCode ierr;
1472bf7f4e0aSPeter Brune   PetscFunctionBegin;
1473bf7f4e0aSPeter Brune   if (linesearch->vec_sol) {
1474bf7f4e0aSPeter Brune     ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr);
1475bf7f4e0aSPeter Brune   } else {
1476bf7f4e0aSPeter Brune     SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
1477bf7f4e0aSPeter Brune   }
1478bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1479bf7f4e0aSPeter Brune }
1480bf7f4e0aSPeter Brune 
14816a388c36SPeter Brune 
1482bf7f4e0aSPeter Brune #undef __FUNCT__
1483f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetSuccess"
1484f40b411bSPeter Brune /*@
1485f1c6b773SPeter Brune    SNESLineSearchGetSuccess - Gets the success/failure status of the last line search application
1486f40b411bSPeter Brune 
1487f40b411bSPeter Brune    Input Parameters:
148878bcb3b5SPeter Brune .  linesearch - linesearch context
1489f40b411bSPeter Brune 
1490f40b411bSPeter Brune    Output Parameters:
14918e557f58SPeter Brune .  success - The success or failure status
1492f40b411bSPeter Brune 
1493cd7522eaSPeter Brune    Notes:
1494cd7522eaSPeter Brune    This is typically called after SNESLineSearchApply in order to determine if the line-search failed
1495cd7522eaSPeter Brune    (and set the SNES convergence accordingly).
1496cd7522eaSPeter Brune 
1497f40b411bSPeter Brune    Level: intermediate
1498f40b411bSPeter Brune 
1499f1c6b773SPeter Brune .seealso: SNESLineSearchSetSuccess()
1500f40b411bSPeter Brune @*/
1501f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetSuccess(SNESLineSearch linesearch, PetscBool *success)
1502bf7f4e0aSPeter Brune {
1503bf7f4e0aSPeter Brune   PetscFunctionBegin;
1504f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
15056a388c36SPeter Brune   PetscValidPointer(success, 2);
1506bf7f4e0aSPeter Brune   if (success) {
1507bf7f4e0aSPeter Brune     *success = linesearch->success;
1508bf7f4e0aSPeter Brune   }
1509bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1510bf7f4e0aSPeter Brune }
1511bf7f4e0aSPeter Brune 
1512bf7f4e0aSPeter Brune #undef __FUNCT__
1513f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetSuccess"
1514f40b411bSPeter Brune /*@
1515f1c6b773SPeter Brune    SNESLineSearchSetSuccess - Sets the success/failure status of the last line search application
1516f40b411bSPeter Brune 
1517f40b411bSPeter Brune    Input Parameters:
151878bcb3b5SPeter Brune +  linesearch - linesearch context
15198e557f58SPeter Brune -  success - The success or failure status
1520f40b411bSPeter Brune 
1521cd7522eaSPeter Brune    Notes:
1522cd7522eaSPeter Brune    This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set
1523cd7522eaSPeter Brune    the success or failure of the line search method.
1524cd7522eaSPeter Brune 
152578bcb3b5SPeter Brune    Level: developer
1526f40b411bSPeter Brune 
1527f1c6b773SPeter Brune .seealso: SNESLineSearchGetSuccess()
1528f40b411bSPeter Brune @*/
1529f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetSuccess(SNESLineSearch linesearch, PetscBool success)
15306a388c36SPeter Brune {
1531f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
15326a388c36SPeter Brune   PetscFunctionBegin;
15336a388c36SPeter Brune   linesearch->success = success;
15346a388c36SPeter Brune   PetscFunctionReturn(0);
15356a388c36SPeter Brune }
15366a388c36SPeter Brune 
15376a388c36SPeter Brune #undef __FUNCT__
1538f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetVIFunctions"
15399bd66eb0SPeter Brune /*@C
1540f1c6b773SPeter Brune    SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.
15419bd66eb0SPeter Brune 
15429bd66eb0SPeter Brune    Input Parameters:
15439bd66eb0SPeter Brune +  snes - nonlinear context obtained from SNESCreate()
15449bd66eb0SPeter Brune .  projectfunc - function for projecting the function to the bounds
15459bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
15469bd66eb0SPeter Brune 
15479bd66eb0SPeter Brune    Logically Collective on SNES
15489bd66eb0SPeter Brune 
15499bd66eb0SPeter Brune    Calling sequence of projectfunc:
15509bd66eb0SPeter Brune .vb
15519bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X)
15529bd66eb0SPeter Brune .ve
15539bd66eb0SPeter Brune 
15549bd66eb0SPeter Brune     Input parameters for projectfunc:
15559bd66eb0SPeter Brune +   snes - nonlinear context
15569bd66eb0SPeter Brune -   X - current solution
15579bd66eb0SPeter Brune 
1558cd7522eaSPeter Brune     Output parameters for projectfunc:
15599bd66eb0SPeter Brune .   X - Projected solution
15609bd66eb0SPeter Brune 
15619bd66eb0SPeter Brune    Calling sequence of normfunc:
15629bd66eb0SPeter Brune .vb
15639bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm)
15649bd66eb0SPeter Brune .ve
15659bd66eb0SPeter Brune 
1566cd7522eaSPeter Brune     Input parameters for normfunc:
15679bd66eb0SPeter Brune +   snes - nonlinear context
15689bd66eb0SPeter Brune .   X - current solution
15699bd66eb0SPeter Brune -   F - current residual
15709bd66eb0SPeter Brune 
1571cd7522eaSPeter Brune     Output parameters for normfunc:
15729bd66eb0SPeter Brune .   fnorm - VI-specific norm of the function
15739bd66eb0SPeter Brune 
1574cd7522eaSPeter Brune     Notes:
1575cd7522eaSPeter Brune     The VI solvers require projection of the solution to the feasible set.  projectfunc should implement this.
1576cd7522eaSPeter Brune 
1577cd7522eaSPeter Brune     The VI solvers require special evaluation of the function norm such that the norm is only calculated
1578cd7522eaSPeter Brune     on the inactive set.  This should be implemented by normfunc.
15799bd66eb0SPeter Brune 
15809bd66eb0SPeter Brune     Level: developer
15819bd66eb0SPeter Brune 
15829bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, set, line search
15839bd66eb0SPeter Brune 
1584f1c6b773SPeter Brune .seealso: SNESLineSearchGetVIFunctions(), SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck()
15859bd66eb0SPeter Brune @*/
1586f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc)
15879bd66eb0SPeter Brune {
15889bd66eb0SPeter Brune   PetscFunctionBegin;
1589f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
15909bd66eb0SPeter Brune   if (projectfunc) linesearch->ops->viproject = projectfunc;
15919bd66eb0SPeter Brune   if (normfunc) linesearch->ops->vinorm = normfunc;
15929bd66eb0SPeter Brune   PetscFunctionReturn(0);
15939bd66eb0SPeter Brune }
15949bd66eb0SPeter Brune 
15959bd66eb0SPeter Brune #undef __FUNCT__
1596f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetVIFunctions"
15979bd66eb0SPeter Brune /*@C
1598f1c6b773SPeter Brune    SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.
15999bd66eb0SPeter Brune 
16009bd66eb0SPeter Brune    Input Parameters:
16019bd66eb0SPeter Brune .  snes - nonlinear context obtained from SNESCreate()
16029bd66eb0SPeter Brune 
16039bd66eb0SPeter Brune    Output Parameters:
16049bd66eb0SPeter Brune +  projectfunc - function for projecting the function to the bounds
16059bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
16069bd66eb0SPeter Brune 
16079bd66eb0SPeter Brune    Logically Collective on SNES
16089bd66eb0SPeter Brune 
16099bd66eb0SPeter Brune     Level: developer
16109bd66eb0SPeter Brune 
16119bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, get, line search
16129bd66eb0SPeter Brune 
1613f1c6b773SPeter Brune .seealso: SNESLineSearchSetVIFunctions(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
16149bd66eb0SPeter Brune @*/
1615f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc)
16169bd66eb0SPeter Brune {
16179bd66eb0SPeter Brune   PetscFunctionBegin;
16189bd66eb0SPeter Brune   if (projectfunc) *projectfunc = linesearch->ops->viproject;
16199bd66eb0SPeter Brune   if (normfunc) *normfunc = linesearch->ops->vinorm;
16209bd66eb0SPeter Brune   PetscFunctionReturn(0);
16219bd66eb0SPeter Brune }
16229bd66eb0SPeter Brune 
16239bd66eb0SPeter Brune #undef __FUNCT__
1624f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchRegister"
1625bf7f4e0aSPeter Brune /*@C
1626f1c6b773SPeter Brune   SNESLineSearchRegister - See SNESLineSearchRegisterDynamic()
1627bf7f4e0aSPeter Brune 
1628bf7f4e0aSPeter Brune   Level: advanced
1629bf7f4e0aSPeter Brune @*/
1630f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNESLineSearch))
1631bf7f4e0aSPeter Brune {
1632bf7f4e0aSPeter Brune   char           fullname[PETSC_MAX_PATH_LEN];
1633bf7f4e0aSPeter Brune   PetscErrorCode ierr;
1634bf7f4e0aSPeter Brune 
1635bf7f4e0aSPeter Brune   PetscFunctionBegin;
1636bf7f4e0aSPeter Brune   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
1637f1c6b773SPeter Brune   ierr = PetscFListAdd(&SNESLineSearchList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
1638bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1639bf7f4e0aSPeter Brune }
1640