xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision 3c7d666354b1ae9ca8c08ab255ad3ef184f4ea79)
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,
36*3c7d6663SPeter Brune                            "SNESLineSearch","Linesearch","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:
2907b1df9c1SPeter Brune +  linesearch - The linesearch instance.
2917b1df9c1SPeter Brune .  X - The current solution
2927b1df9c1SPeter Brune -  Y - The step direction
293f40b411bSPeter Brune 
294f40b411bSPeter Brune    Output Parameters:
2958e557f58SPeter Brune .  changed - Indicator that the precheck routine has changed anything
296f40b411bSPeter Brune 
297f40b411bSPeter Brune    Level: Beginner
298f40b411bSPeter Brune 
299f1c6b773SPeter Brune .keywords: SNESLineSearch, Create
300f40b411bSPeter Brune 
301f1c6b773SPeter Brune .seealso: SNESLineSearchPostCheck()
302f40b411bSPeter Brune @*/
3037b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed)
304bf7f4e0aSPeter Brune {
305bf7f4e0aSPeter Brune   PetscErrorCode ierr;
306bf7f4e0aSPeter Brune   PetscFunctionBegin;
307bf7f4e0aSPeter Brune   *changed = PETSC_FALSE;
308bf7f4e0aSPeter Brune   if (linesearch->ops->precheckstep) {
3097b1df9c1SPeter Brune     ierr = (*linesearch->ops->precheckstep)(linesearch, X, Y, changed, linesearch->precheckctx);CHKERRQ(ierr);
310bf7f4e0aSPeter Brune   }
311bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
312bf7f4e0aSPeter Brune }
313bf7f4e0aSPeter Brune 
314bf7f4e0aSPeter Brune #undef __FUNCT__
315f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPostCheck"
316f40b411bSPeter Brune /*@
317f1c6b773SPeter Brune    SNESLineSearchPostCheck - Prepares the line search for being applied.
318f40b411bSPeter Brune 
319cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
320f40b411bSPeter Brune 
321f40b411bSPeter Brune    Input Parameters:
3227b1df9c1SPeter Brune +  linesearch - The linesearch context
3237b1df9c1SPeter Brune .  X - The last solution
3247b1df9c1SPeter Brune .  Y - The step direction
3257b1df9c1SPeter Brune -  W - The updated solution, W = X + lambda*Y for some lambda
326f40b411bSPeter Brune 
327f40b411bSPeter Brune    Output Parameters:
32878bcb3b5SPeter Brune +  changed_Y - Indicator if the direction Y has been changed.
32978bcb3b5SPeter Brune -  changed_W - Indicator if the new candidate solution W has been changed.
330f40b411bSPeter Brune 
331f40b411bSPeter Brune    Level: Intermediate
332f40b411bSPeter Brune 
333f1c6b773SPeter Brune .keywords: SNESLineSearch, Create
334f40b411bSPeter Brune 
335f1c6b773SPeter Brune .seealso: SNESLineSearchPreCheck()
336f40b411bSPeter Brune @*/
3377b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch,Vec X,Vec Y,Vec W,PetscBool *changed_Y,PetscBool *changed_W)
338bf7f4e0aSPeter Brune {
339bf7f4e0aSPeter Brune   PetscErrorCode ierr;
340bf7f4e0aSPeter Brune   PetscFunctionBegin;
341bf7f4e0aSPeter Brune   *changed_Y = PETSC_FALSE;
342bf7f4e0aSPeter Brune   *changed_W = PETSC_FALSE;
343bf7f4e0aSPeter Brune   if (linesearch->ops->postcheckstep) {
3447b1df9c1SPeter Brune     ierr = (*linesearch->ops->postcheckstep)(linesearch,X,Y,W,changed_Y,changed_W,linesearch->postcheckctx);CHKERRQ(ierr);
34586d74e61SPeter Brune   }
34686d74e61SPeter Brune   PetscFunctionReturn(0);
34786d74e61SPeter Brune }
34886d74e61SPeter Brune 
34986d74e61SPeter Brune 
35086d74e61SPeter Brune #undef __FUNCT__
351f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPreCheckPicard"
35286d74e61SPeter Brune /*@C
35386d74e61SPeter Brune    SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration
35486d74e61SPeter Brune 
355cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
35686d74e61SPeter Brune 
35786d74e61SPeter Brune    Input Arguments:
35886d74e61SPeter Brune +  linesearch - linesearch context
35986d74e61SPeter Brune .  X - base state for this step
36086d74e61SPeter Brune .  Y - initial correction
36186d74e61SPeter Brune 
36286d74e61SPeter Brune    Output Arguments:
36386d74e61SPeter Brune +  Y - correction, possibly modified
36486d74e61SPeter Brune -  changed - flag indicating that Y was modified
36586d74e61SPeter Brune 
36686d74e61SPeter Brune    Options Database Key:
367cd7522eaSPeter Brune +  -snes_linesearch_precheck_picard - activate this routine
368cd7522eaSPeter Brune -  -snes_linesearch_precheck_picard_angle - angle
36986d74e61SPeter Brune 
37086d74e61SPeter Brune    Level: advanced
37186d74e61SPeter Brune 
37286d74e61SPeter Brune    Notes:
37386d74e61SPeter Brune    This function should be passed to SNESLineSearchSetPreCheck()
37486d74e61SPeter Brune 
37586d74e61SPeter Brune    The justification for this method involves the linear convergence of a Picard iteration
37686d74e61SPeter Brune    so the Picard linearization should be provided in place of the "Jacobian". This correction
37786d74e61SPeter Brune    is generally not useful when using a Newton linearization.
37886d74e61SPeter Brune 
37986d74e61SPeter Brune    Reference:
38086d74e61SPeter Brune    Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology.
38186d74e61SPeter Brune 
38286d74e61SPeter Brune .seealso: SNESLineSearchSetPreCheck()
38386d74e61SPeter Brune @*/
384f1c6b773SPeter Brune PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed,void *ctx)
38586d74e61SPeter Brune {
38686d74e61SPeter Brune   PetscErrorCode ierr;
38786d74e61SPeter Brune   PetscReal      angle = *(PetscReal*)linesearch->precheckctx;
38886d74e61SPeter Brune   Vec            Ylast;
38986d74e61SPeter Brune   PetscScalar    dot;
390c87759e9SPeter Brune 
39186d74e61SPeter Brune   PetscInt       iter;
39286d74e61SPeter Brune   PetscReal      ynorm,ylastnorm,theta,angle_radians;
39386d74e61SPeter Brune   SNES           snes;
39486d74e61SPeter Brune 
39586d74e61SPeter Brune   PetscFunctionBegin;
396f1c6b773SPeter Brune   ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
39786d74e61SPeter Brune   ierr = PetscObjectQuery((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject*)&Ylast);CHKERRQ(ierr);
39886d74e61SPeter Brune   if (!Ylast) {
39986d74e61SPeter Brune     ierr = VecDuplicate(Y,&Ylast);CHKERRQ(ierr);
40086d74e61SPeter Brune     ierr = PetscObjectCompose((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject)Ylast);CHKERRQ(ierr);
40186d74e61SPeter Brune     ierr = PetscObjectDereference((PetscObject)Ylast);CHKERRQ(ierr);
40286d74e61SPeter Brune   }
40386d74e61SPeter Brune   ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr);
40486d74e61SPeter Brune   if (iter < 2) {
40586d74e61SPeter Brune     ierr = VecCopy(Y,Ylast);CHKERRQ(ierr);
40686d74e61SPeter Brune     *changed = PETSC_FALSE;
40786d74e61SPeter Brune     PetscFunctionReturn(0);
40886d74e61SPeter Brune   }
40986d74e61SPeter Brune 
41086d74e61SPeter Brune   ierr = VecDot(Y,Ylast,&dot);CHKERRQ(ierr);
41186d74e61SPeter Brune   ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr);
41286d74e61SPeter Brune   ierr = VecNorm(Ylast,NORM_2,&ylastnorm);CHKERRQ(ierr);
41386d74e61SPeter Brune   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
41486d74e61SPeter Brune   theta = acos((double)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0));
41586d74e61SPeter Brune   angle_radians = angle * PETSC_PI / 180.;
41686d74e61SPeter Brune   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
41786d74e61SPeter Brune     /* Modify the step Y */
41886d74e61SPeter Brune     PetscReal alpha,ydiffnorm;
41986d74e61SPeter Brune     ierr = VecAXPY(Ylast,-1.0,Y);CHKERRQ(ierr);
42086d74e61SPeter Brune     ierr = VecNorm(Ylast,NORM_2,&ydiffnorm);CHKERRQ(ierr);
42186d74e61SPeter Brune     alpha = ylastnorm / ydiffnorm;
42286d74e61SPeter Brune     ierr = VecCopy(Y,Ylast);CHKERRQ(ierr);
42386d74e61SPeter Brune     ierr = VecScale(Y,alpha);CHKERRQ(ierr);
42486d74e61SPeter 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);
42586d74e61SPeter Brune   } else {
42686d74e61SPeter Brune     ierr = PetscInfo2(snes,"Angle %G degrees exceeds threshold %G, no correction applied\n",theta*180./PETSC_PI,angle);CHKERRQ(ierr);
42786d74e61SPeter Brune     ierr = VecCopy(Y,Ylast);CHKERRQ(ierr);
42886d74e61SPeter Brune     *changed = PETSC_FALSE;
429bf7f4e0aSPeter Brune   }
430bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
431bf7f4e0aSPeter Brune }
432bf7f4e0aSPeter Brune 
433bf7f4e0aSPeter Brune #undef __FUNCT__
434f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchApply"
435f40b411bSPeter Brune /*@
436cd7522eaSPeter Brune    SNESLineSearchApply - Computes the line-search update.
437f40b411bSPeter Brune 
438f1c6b773SPeter Brune    Collective on SNESLineSearch
439f40b411bSPeter Brune 
440f40b411bSPeter Brune    Input Parameters:
4418e557f58SPeter Brune +  linesearch - The linesearch context
4428e557f58SPeter Brune .  X - The current solution
4438e557f58SPeter Brune .  F - The current function
4448e557f58SPeter Brune .  fnorm - The current norm
4458e557f58SPeter Brune -  Y - The search direction
446f40b411bSPeter Brune 
447f40b411bSPeter Brune    Output Parameters:
4488e557f58SPeter Brune +  X - The new solution
4498e557f58SPeter Brune .  F - The new function
4508e557f58SPeter Brune -  fnorm - The new function norm
451f40b411bSPeter Brune 
452cd7522eaSPeter Brune    Options Database Keys:
453*3c7d6663SPeter Brune + -snes_linesearch_type - basic, bt, l2, cp, shell
454cd7522eaSPeter Brune . -snes_linesearch_monitor - Print progress of line searches
455*3c7d6663SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter
456cd7522eaSPeter Brune . -snes_linesearch_norms   - Turn on/off the linesearch norms
457*3c7d6663SPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess
458*3c7d6663SPeter Brune - -snes_linesearch_max_it - The number of iterations for iterative line searches
459cd7522eaSPeter Brune 
460cd7522eaSPeter Brune    Notes:
461cd7522eaSPeter Brune    This is typically called from within a SNESSolve() implementation in order to
462cd7522eaSPeter Brune    help with convergence of the nonlinear method.  Various SNES types use line searches
463cd7522eaSPeter Brune    in different ways, but the overarching theme is that a line search is used to determine
464cd7522eaSPeter Brune    an optimal damping parameter of a step at each iteration of the method.  Each
465cd7522eaSPeter Brune    application of the line search may invoke SNESComputeFunction several times, and
466cd7522eaSPeter Brune    therefore may be fairly expensive.
467cd7522eaSPeter Brune 
468f40b411bSPeter Brune    Level: Intermediate
469f40b411bSPeter Brune 
470f1c6b773SPeter Brune .keywords: SNESLineSearch, Create
471f40b411bSPeter Brune 
472cd7522eaSPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchPreCheck(), SNESLineSearchPostCheck(), SNESSolve(), SNESComputeFunction()
473f40b411bSPeter Brune @*/
474f1c6b773SPeter Brune PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y) {
475bf7f4e0aSPeter Brune   PetscErrorCode ierr;
476bf7f4e0aSPeter Brune   PetscFunctionBegin;
477bf7f4e0aSPeter Brune 
478bf7f4e0aSPeter Brune   /* check the pointers */
479f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
480bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
481bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
482bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
483bf7f4e0aSPeter Brune 
484bf7f4e0aSPeter Brune   linesearch->success = PETSC_TRUE;
485bf7f4e0aSPeter Brune 
486bf7f4e0aSPeter Brune   linesearch->vec_sol = X;
487bf7f4e0aSPeter Brune   linesearch->vec_update = Y;
488bf7f4e0aSPeter Brune   linesearch->vec_func = F;
489bf7f4e0aSPeter Brune 
490f1c6b773SPeter Brune   ierr = SNESLineSearchSetUp(linesearch);CHKERRQ(ierr);
491bf7f4e0aSPeter Brune 
492bf7f4e0aSPeter Brune   if (!linesearch->keeplambda)
493bf7f4e0aSPeter Brune     linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */
494bf7f4e0aSPeter Brune 
495bf7f4e0aSPeter Brune   if (fnorm) {
496bf7f4e0aSPeter Brune     linesearch->fnorm = *fnorm;
497bf7f4e0aSPeter Brune   } else {
498bf7f4e0aSPeter Brune     ierr = VecNorm(F, NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
499bf7f4e0aSPeter Brune   }
500bf7f4e0aSPeter Brune 
501f1c6b773SPeter Brune   ierr = PetscLogEventBegin(SNESLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
502bf7f4e0aSPeter Brune 
503bf7f4e0aSPeter Brune   ierr = (*linesearch->ops->apply)(linesearch);CHKERRQ(ierr);
504bf7f4e0aSPeter Brune 
505f1c6b773SPeter Brune   ierr = PetscLogEventEnd(SNESLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
506bf7f4e0aSPeter Brune 
507bf7f4e0aSPeter Brune   if (fnorm)
508bf7f4e0aSPeter Brune     *fnorm = linesearch->fnorm;
509bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
510bf7f4e0aSPeter Brune }
511bf7f4e0aSPeter Brune 
512bf7f4e0aSPeter Brune #undef __FUNCT__
513f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchDestroy"
514f40b411bSPeter Brune /*@
515f1c6b773SPeter Brune    SNESLineSearchDestroy - Destroys the line search instance.
516f40b411bSPeter Brune 
517f1c6b773SPeter Brune    Collective on SNESLineSearch
518f40b411bSPeter Brune 
519f40b411bSPeter Brune    Input Parameters:
5208e557f58SPeter Brune .  linesearch - The linesearch context
521f40b411bSPeter Brune 
522f40b411bSPeter Brune    Level: Intermediate
523f40b411bSPeter Brune 
52478bcb3b5SPeter Brune .keywords: SNESLineSearch, Destroy
525f40b411bSPeter Brune 
526cd7522eaSPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchReset(), SNESDestroy()
527f40b411bSPeter Brune @*/
528f1c6b773SPeter Brune PetscErrorCode SNESLineSearchDestroy(SNESLineSearch * linesearch) {
529bf7f4e0aSPeter Brune   PetscErrorCode ierr;
530bf7f4e0aSPeter Brune   PetscFunctionBegin;
531bf7f4e0aSPeter Brune   if (!*linesearch) PetscFunctionReturn(0);
532f1c6b773SPeter Brune   PetscValidHeaderSpecific((*linesearch),SNESLINESEARCH_CLASSID,1);
533bf7f4e0aSPeter Brune   if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = 0; PetscFunctionReturn(0);}
534bf7f4e0aSPeter Brune   ierr = PetscObjectDepublish((*linesearch));CHKERRQ(ierr);
535f1c6b773SPeter Brune   ierr = SNESLineSearchReset(*linesearch);
536bf7f4e0aSPeter Brune   if ((*linesearch)->ops->destroy) {
537bf7f4e0aSPeter Brune     (*linesearch)->ops->destroy(*linesearch);
538bf7f4e0aSPeter Brune   }
539bf7f4e0aSPeter Brune   ierr = PetscViewerDestroy(&(*linesearch)->monitor);CHKERRQ(ierr);
540e7058c64SPeter Brune   ierr = PetscHeaderDestroy(linesearch);CHKERRQ(ierr);
541bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
542bf7f4e0aSPeter Brune }
543bf7f4e0aSPeter Brune 
544bf7f4e0aSPeter Brune #undef __FUNCT__
545f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetMonitor"
546f40b411bSPeter Brune /*@
547cd7522eaSPeter Brune    SNESLineSearchSetMonitor - Turns on/off printing useful information and debugging output about the line search.
548bf7f4e0aSPeter Brune 
549bf7f4e0aSPeter Brune    Input Parameters:
550bf7f4e0aSPeter Brune +  snes - nonlinear context obtained from SNESCreate()
551bf7f4e0aSPeter Brune -  flg - PETSC_TRUE to monitor the line search
552bf7f4e0aSPeter Brune 
553bf7f4e0aSPeter Brune    Logically Collective on SNES
554bf7f4e0aSPeter Brune 
555bf7f4e0aSPeter Brune    Options Database:
5568e557f58SPeter Brune .   -snes_linesearch_monitor - enables the monitor
557bf7f4e0aSPeter Brune 
558bf7f4e0aSPeter Brune    Level: intermediate
559bf7f4e0aSPeter Brune 
560bf7f4e0aSPeter Brune 
561cd7522eaSPeter Brune .seealso: SNESLineSearchGetMonitor(), PetscViewer
562bf7f4e0aSPeter Brune @*/
563f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetMonitor(SNESLineSearch linesearch, PetscBool flg)
564bf7f4e0aSPeter Brune {
565bf7f4e0aSPeter Brune 
566bf7f4e0aSPeter Brune   PetscErrorCode ierr;
567bf7f4e0aSPeter Brune   PetscFunctionBegin;
568bf7f4e0aSPeter Brune   if (flg && !linesearch->monitor) {
569bf7f4e0aSPeter Brune     ierr = PetscViewerASCIIOpen(((PetscObject)linesearch)->comm,"stdout",&linesearch->monitor);CHKERRQ(ierr);
570bf7f4e0aSPeter Brune   } else if (!flg && linesearch->monitor) {
571bf7f4e0aSPeter Brune     ierr = PetscViewerDestroy(&linesearch->monitor);CHKERRQ(ierr);
572bf7f4e0aSPeter Brune   }
573bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
574bf7f4e0aSPeter Brune }
575bf7f4e0aSPeter Brune 
576bf7f4e0aSPeter Brune #undef __FUNCT__
577f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetMonitor"
578f40b411bSPeter Brune /*@
579cd7522eaSPeter Brune    SNESLineSearchGetMonitor - Gets the PetscViewer instance for the line search monitor.
5806a388c36SPeter Brune 
581f40b411bSPeter Brune    Input Parameters:
5828e557f58SPeter Brune .  linesearch - linesearch context
583f40b411bSPeter Brune 
584f40b411bSPeter Brune    Input Parameters:
5858e557f58SPeter Brune .  monitor - monitor context
586f40b411bSPeter Brune 
587f40b411bSPeter Brune    Logically Collective on SNES
588f40b411bSPeter Brune 
589f40b411bSPeter Brune 
590f40b411bSPeter Brune    Options Database Keys:
5918e557f58SPeter Brune .   -snes_linesearch_monitor - enables the monitor
592f40b411bSPeter Brune 
593f40b411bSPeter Brune    Level: intermediate
594f40b411bSPeter Brune 
595f40b411bSPeter Brune 
596cd7522eaSPeter Brune .seealso: SNESLineSearchSetMonitor(), PetscViewer
597f40b411bSPeter Brune @*/
598f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetMonitor(SNESLineSearch linesearch, PetscViewer *monitor)
5996a388c36SPeter Brune {
6006a388c36SPeter Brune 
6016a388c36SPeter Brune   PetscFunctionBegin;
602f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
6036a388c36SPeter Brune   if (monitor) {
6046a388c36SPeter Brune     PetscValidPointer(monitor, 2);
6056a388c36SPeter Brune     *monitor = linesearch->monitor;
6066a388c36SPeter Brune   }
6076a388c36SPeter Brune   PetscFunctionReturn(0);
6086a388c36SPeter Brune }
6096a388c36SPeter Brune 
6106a388c36SPeter Brune #undef __FUNCT__
611f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetFromOptions"
612f40b411bSPeter Brune /*@
613f1c6b773SPeter Brune    SNESLineSearchSetFromOptions - Sets options for the line search
614f40b411bSPeter Brune 
615f40b411bSPeter Brune    Input Parameters:
6168e557f58SPeter Brune .  linesearch - linesearch context
617f40b411bSPeter Brune 
618f40b411bSPeter Brune    Options Database Keys:
619*3c7d6663SPeter Brune + -snes_linesearch_type <type> - basic, bt, l2, cp, shell
620*3c7d6663SPeter Brune . -snes_linesearch_order <order> - 1, 2, 3.  Most types only support certain orders (bt supports 2 or 3)
62171eef1aeSPeter Brune . -snes_linesearch_minlambda - The minimum step length
6221a9b3a06SPeter Brune . -snes_linesearch_maxstep - The maximum step size
6231a9b3a06SPeter Brune . -snes_linesearch_rtol - Relative tolerance for iterative line searches
6241a9b3a06SPeter Brune . -snes_linesearch_atol - Absolute tolerance for iterative line searches
6251a9b3a06SPeter Brune . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches
6261a9b3a06SPeter Brune . -snes_linesearch_max_it - The number of iterations for iterative line searches
627cd7522eaSPeter Brune . -snes_linesearch_monitor - Print progress of line searches
6288e557f58SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter
629cd7522eaSPeter Brune . -snes_linesearch_norms   - Turn on/off the linesearch norms
630cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess.
6311a9b3a06SPeter Brune . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method
6321a9b3a06SPeter Brune - -snes_linesearch_precheck_picard_angle - Angle used in picard precheck method
633f40b411bSPeter Brune 
634f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
635f40b411bSPeter Brune 
636f40b411bSPeter Brune    Level: intermediate
637f40b411bSPeter Brune 
638f40b411bSPeter Brune 
639*3c7d6663SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetOrder(), SNESLineSearchSetType(), SNESLineSearchSetTolerances(), SNESLineSearchSetDamping(), SNESLineSearchPreCheckPicard()
640f40b411bSPeter Brune @*/
641f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch) {
642bf7f4e0aSPeter Brune   PetscErrorCode ierr;
6431a4f838cSPeter Brune   const char     *deft = SNESLINESEARCHBASIC;
644bf7f4e0aSPeter Brune   char           type[256];
645bf7f4e0aSPeter Brune   PetscBool      flg, set;
646bf7f4e0aSPeter Brune   PetscFunctionBegin;
647f1c6b773SPeter Brune   if (!SNESLineSearchRegisterAllCalled) {ierr = SNESLineSearchRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
648bf7f4e0aSPeter Brune 
649bf7f4e0aSPeter Brune   ierr = PetscObjectOptionsBegin((PetscObject)linesearch);CHKERRQ(ierr);
650bf7f4e0aSPeter Brune   if (((PetscObject)linesearch)->type_name) {
651bf7f4e0aSPeter Brune     deft = ((PetscObject)linesearch)->type_name;
652bf7f4e0aSPeter Brune   }
653*3c7d6663SPeter Brune   ierr = PetscOptionsList("-snes_linesearch_type","Linesearch type","SNESLineSearchSetType",SNESLineSearchList,deft,type,256,&flg);CHKERRQ(ierr);
654bf7f4e0aSPeter Brune   if (flg) {
655f1c6b773SPeter Brune     ierr = SNESLineSearchSetType(linesearch,type);CHKERRQ(ierr);
656bf7f4e0aSPeter Brune   } else if (!((PetscObject)linesearch)->type_name) {
657f1c6b773SPeter Brune     ierr = SNESLineSearchSetType(linesearch,deft);CHKERRQ(ierr);
658bf7f4e0aSPeter Brune   }
659bf7f4e0aSPeter Brune   if (linesearch->ops->setfromoptions) {
660bf7f4e0aSPeter Brune     (*linesearch->ops->setfromoptions)(linesearch);CHKERRQ(ierr);
661bf7f4e0aSPeter Brune   }
662bf7f4e0aSPeter Brune 
6637a35526eSPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_monitor","Print progress of line searches","SNESSNESLineSearchSetMonitor",
664bf7f4e0aSPeter Brune                           linesearch->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
665f1c6b773SPeter Brune   if (set) {ierr = SNESLineSearchSetMonitor(linesearch,flg);CHKERRQ(ierr);}
666bf7f4e0aSPeter Brune 
6671a9b3a06SPeter Brune   /* tolerances */
66871eef1aeSPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_minlambda","Minimum step length","SNESLineSearchSetTolerances",linesearch->steptol,&linesearch->steptol,0);CHKERRQ(ierr);
6691a9b3a06SPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_maxstep","Maximum step size","SNESLineSearchSetTolerances",linesearch->maxstep,&linesearch->maxstep,0);CHKERRQ(ierr);
6701a9b3a06SPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_rtol","Relative tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->rtol,&linesearch->rtol,0);CHKERRQ(ierr);
6711a9b3a06SPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_atol","Absolute tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->atol,&linesearch->atol,0);CHKERRQ(ierr);
6721a9b3a06SPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_ltol","Change in lambda tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->ltol,&linesearch->ltol,0);CHKERRQ(ierr);
6738e557f58SPeter Brune   ierr = PetscOptionsInt("-snes_linesearch_max_it","Maximum iterations for iterative line searches","SNESLineSearchSetTolerances",linesearch->max_its,&linesearch->max_its,0);CHKERRQ(ierr);
6747a35526eSPeter Brune 
6751a9b3a06SPeter Brune   /* damping parameters */
6761a9b3a06SPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_damping","Line search damping and initial step guess","SNESLineSearchSetDamping",linesearch->damping,&linesearch->damping,0);CHKERRQ(ierr);
6771a9b3a06SPeter Brune 
6781a9b3a06SPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_keeplambda","Use previous lambda as damping","SNESLineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,0);CHKERRQ(ierr);
6791a9b3a06SPeter Brune 
6801a9b3a06SPeter Brune   /* precheck */
6817a35526eSPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_precheck_picard","Use a correction that sometimes improves convergence of Picard iteration","SNESLineSearchPreCheckPicard",flg,&flg,&set);CHKERRQ(ierr);
6827a35526eSPeter Brune   if (set) {
6837a35526eSPeter Brune     if (flg) {
6847a35526eSPeter Brune       linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */
6857a35526eSPeter Brune       ierr = PetscOptionsReal("-snes_linesearch_precheck_picard_angle","Maximum angle at which to activate the correction",
6867a35526eSPeter Brune                               "none",linesearch->precheck_picard_angle,&linesearch->precheck_picard_angle,PETSC_NULL);CHKERRQ(ierr);
6877a35526eSPeter Brune       ierr = SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&linesearch->precheck_picard_angle);CHKERRQ(ierr);
6887a35526eSPeter Brune     } else {
6897a35526eSPeter Brune       ierr = SNESLineSearchSetPreCheck(linesearch,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
6907a35526eSPeter Brune     }
6917a35526eSPeter Brune   }
692b000cd8dSPeter Brune   ierr = PetscOptionsInt("-snes_linesearch_order","Order of approximation used in the line search","SNESLineSearchSetOrder",linesearch->order,&linesearch->order,0);CHKERRQ(ierr);
6931a9b3a06SPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_norms","Compute final norms in line search","SNESLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,0);CHKERRQ(ierr);
6947a35526eSPeter Brune 
695bf7f4e0aSPeter Brune   ierr = PetscObjectProcessOptionsHandlers((PetscObject)linesearch);CHKERRQ(ierr);
696bf7f4e0aSPeter Brune   ierr = PetscOptionsEnd();CHKERRQ(ierr);
697bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
698bf7f4e0aSPeter Brune }
699bf7f4e0aSPeter Brune 
700bf7f4e0aSPeter Brune #undef __FUNCT__
701f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchView"
702f40b411bSPeter Brune /*@
703cd7522eaSPeter Brune    SNESLineSearchView - Prints useful information about the line search not
704cd7522eaSPeter Brune    related to an individual call.
705f40b411bSPeter Brune 
706f40b411bSPeter Brune    Input Parameters:
7078e557f58SPeter Brune .  linesearch - linesearch context
708f40b411bSPeter Brune 
709f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
710f40b411bSPeter Brune 
711f40b411bSPeter Brune    Level: intermediate
712f40b411bSPeter Brune 
713f1c6b773SPeter Brune .seealso: SNESLineSearchCreate()
714f40b411bSPeter Brune @*/
7157f1410a3SPeter Brune PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer) {
7167f1410a3SPeter Brune   PetscErrorCode ierr;
7177f1410a3SPeter Brune   PetscBool      iascii;
718bf7f4e0aSPeter Brune   PetscFunctionBegin;
7197f1410a3SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
7207f1410a3SPeter Brune   if (!viewer) {
7217f1410a3SPeter Brune     ierr = PetscViewerASCIIGetStdout(((PetscObject)linesearch)->comm,&viewer);CHKERRQ(ierr);
7227f1410a3SPeter Brune   }
7237f1410a3SPeter Brune   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
7247f1410a3SPeter Brune   PetscCheckSameComm(linesearch,1,viewer,2);
725f40b411bSPeter Brune 
726251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
7277f1410a3SPeter Brune   if (iascii) {
7287f1410a3SPeter Brune     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)linesearch,viewer,"SNESLineSearch Object");CHKERRQ(ierr);
7297f1410a3SPeter Brune     if (linesearch->ops->view) {
7307f1410a3SPeter Brune       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
7317f1410a3SPeter Brune       ierr = (*linesearch->ops->view)(linesearch,viewer);CHKERRQ(ierr);
7327f1410a3SPeter Brune       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
7337f1410a3SPeter Brune     }
734007b6e36SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  maxstep=%e, minlambda=%e\n", linesearch->maxstep,linesearch->steptol);CHKERRQ(ierr);
735007b6e36SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%e, absolute=%e, lambda=%e\n", linesearch->rtol,linesearch->atol,linesearch->ltol);CHKERRQ(ierr);
7367f1410a3SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D\n", linesearch->max_its);CHKERRQ(ierr);
7377f1410a3SPeter Brune     if (linesearch->ops->precheckstep) {
7387f1410a3SPeter Brune       if (linesearch->ops->precheckstep == SNESLineSearchPreCheckPicard) {
7397f1410a3SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  using precheck step to speed up Picard convergence\n", linesearch->max_its);CHKERRQ(ierr);
7407f1410a3SPeter Brune       } else {
7417f1410a3SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  using user-defined precheck step\n", linesearch->max_its);CHKERRQ(ierr);
7427f1410a3SPeter Brune       }
7437f1410a3SPeter Brune     }
7447f1410a3SPeter Brune     if (linesearch->ops->postcheckstep) {
7457f1410a3SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer,"  using user-defined postcheck step\n", linesearch->max_its);CHKERRQ(ierr);
7467f1410a3SPeter Brune     }
7477f1410a3SPeter Brune   }
748bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
749bf7f4e0aSPeter Brune }
750bf7f4e0aSPeter Brune 
751bf7f4e0aSPeter Brune #undef __FUNCT__
752f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetType"
753ea5d4fccSPeter Brune /*@C
754f1c6b773SPeter Brune    SNESLineSearchSetType - Sets the linesearch type
755f40b411bSPeter Brune 
756f40b411bSPeter Brune    Input Parameters:
7578e557f58SPeter Brune +  linesearch - linesearch context
758f40b411bSPeter Brune -  type - The type of line search to be used
759f40b411bSPeter Brune 
760cd7522eaSPeter Brune    Available Types:
761cd7522eaSPeter Brune +  basic - Simple damping line search.
7628e557f58SPeter Brune .  bt - Backtracking line search over the L2 norm of the function
7638e557f58SPeter Brune .  l2 - Secant line search over the L2 norm of the function
7648e557f58SPeter Brune .  cp - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x)
7658e557f58SPeter Brune -  shell - User provided SNESLineSearch implementation
766cd7522eaSPeter Brune 
767f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
768f40b411bSPeter Brune 
769f40b411bSPeter Brune    Level: intermediate
770f40b411bSPeter Brune 
771f40b411bSPeter Brune 
772f1c6b773SPeter Brune .seealso: SNESLineSearchCreate()
773f40b411bSPeter Brune @*/
774f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, const SNESLineSearchType type)
775bf7f4e0aSPeter Brune {
776bf7f4e0aSPeter Brune 
777f1c6b773SPeter Brune   PetscErrorCode ierr,(*r)(SNESLineSearch);
778bf7f4e0aSPeter Brune   PetscBool      match;
779bf7f4e0aSPeter Brune 
780bf7f4e0aSPeter Brune   PetscFunctionBegin;
781f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
782bf7f4e0aSPeter Brune   PetscValidCharPointer(type,2);
783bf7f4e0aSPeter Brune 
784251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)linesearch,type,&match);CHKERRQ(ierr);
785bf7f4e0aSPeter Brune   if (match) PetscFunctionReturn(0);
786bf7f4e0aSPeter Brune 
787f1c6b773SPeter Brune   ierr =  PetscFListFind(SNESLineSearchList,((PetscObject)linesearch)->comm,type,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr);
788bf7f4e0aSPeter Brune   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type);
789bf7f4e0aSPeter Brune   /* Destroy the previous private linesearch context */
790bf7f4e0aSPeter Brune   if (linesearch->ops->destroy) {
791bf7f4e0aSPeter Brune     ierr = (*(linesearch)->ops->destroy)(linesearch);CHKERRQ(ierr);
792bf7f4e0aSPeter Brune     linesearch->ops->destroy = PETSC_NULL;
793bf7f4e0aSPeter Brune   }
794f1c6b773SPeter Brune   /* Reinitialize function pointers in SNESLineSearchOps structure */
795bf7f4e0aSPeter Brune   linesearch->ops->apply          = 0;
796bf7f4e0aSPeter Brune   linesearch->ops->view           = 0;
797bf7f4e0aSPeter Brune   linesearch->ops->setfromoptions = 0;
798bf7f4e0aSPeter Brune   linesearch->ops->destroy        = 0;
799bf7f4e0aSPeter Brune 
800bf7f4e0aSPeter Brune   ierr = PetscObjectChangeTypeName((PetscObject)linesearch,type);CHKERRQ(ierr);
801bf7f4e0aSPeter Brune   ierr = (*r)(linesearch);CHKERRQ(ierr);
802bf7f4e0aSPeter Brune #if defined(PETSC_HAVE_AMS)
803bf7f4e0aSPeter Brune   if (PetscAMSPublishAll) {
804bf7f4e0aSPeter Brune     ierr = PetscObjectAMSPublish((PetscObject)linesearch);CHKERRQ(ierr);
805bf7f4e0aSPeter Brune   }
806bf7f4e0aSPeter Brune #endif
807bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
808bf7f4e0aSPeter Brune }
809bf7f4e0aSPeter Brune 
810bf7f4e0aSPeter Brune #undef __FUNCT__
811f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetSNES"
812f40b411bSPeter Brune /*@
81378bcb3b5SPeter Brune    SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation.
814f40b411bSPeter Brune 
815f40b411bSPeter Brune    Input Parameters:
8168e557f58SPeter Brune +  linesearch - linesearch context
817f40b411bSPeter Brune -  snes - The snes instance
818f40b411bSPeter Brune 
81978bcb3b5SPeter Brune    Level: developer
82078bcb3b5SPeter Brune 
82178bcb3b5SPeter Brune    Notes:
82278bcb3b5SPeter Brune    This happens automatically when the line search is gotten/created with
82378bcb3b5SPeter Brune    SNESGetSNESLineSearch().  This routine is therefore mainly called within SNES
82478bcb3b5SPeter Brune    implementations.
825f40b411bSPeter Brune 
8268141a3b9SPeter Brune    Level: developer
827f40b411bSPeter Brune 
828cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
829f40b411bSPeter Brune @*/
830f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes){
831bf7f4e0aSPeter Brune   PetscFunctionBegin;
832f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
833bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,2);
834bf7f4e0aSPeter Brune   linesearch->snes = snes;
835bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
836bf7f4e0aSPeter Brune }
837bf7f4e0aSPeter Brune 
838bf7f4e0aSPeter Brune #undef __FUNCT__
839f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetSNES"
840f40b411bSPeter Brune /*@
8418141a3b9SPeter Brune    SNESLineSearchGetSNES - Gets the SNES instance associated with the line search.
8428141a3b9SPeter Brune    Having an associated SNES is necessary because most line search implementations must be able to
8438141a3b9SPeter Brune    evaluate the function using SNESComputeFunction() for the associated SNES.  This routine
8448141a3b9SPeter Brune    is used in the line search implementations when one must get this associated SNES instance.
845f40b411bSPeter Brune 
846f40b411bSPeter Brune    Input Parameters:
8478e557f58SPeter Brune .  linesearch - linesearch context
848f40b411bSPeter Brune 
849f40b411bSPeter Brune    Output Parameters:
850f40b411bSPeter Brune .  snes - The snes instance
851f40b411bSPeter Brune 
8528141a3b9SPeter Brune    Level: developer
853f40b411bSPeter Brune 
854cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
855f40b411bSPeter Brune @*/
856f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes){
857bf7f4e0aSPeter Brune   PetscFunctionBegin;
858f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
8596a388c36SPeter Brune   PetscValidPointer(snes, 2);
860bf7f4e0aSPeter Brune   *snes = linesearch->snes;
861bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
862bf7f4e0aSPeter Brune }
863bf7f4e0aSPeter Brune 
8646a388c36SPeter Brune #undef __FUNCT__
865f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetLambda"
866f40b411bSPeter Brune /*@
867f1c6b773SPeter Brune    SNESLineSearchGetLambda - Gets the last linesearch steplength discovered.
868f40b411bSPeter Brune 
869f40b411bSPeter Brune    Input Parameters:
8708e557f58SPeter Brune .  linesearch - linesearch context
871f40b411bSPeter Brune 
872f40b411bSPeter Brune    Output Parameters:
873cd7522eaSPeter Brune .  lambda - The last steplength computed during SNESLineSearchApply()
874f40b411bSPeter Brune 
87578bcb3b5SPeter Brune    Level: advanced
87678bcb3b5SPeter Brune 
8778e557f58SPeter Brune    Notes:
8788e557f58SPeter Brune    This is useful in methods where the solver is ill-scaled and
87978bcb3b5SPeter Brune    requires some adaptive notion of the difference in scale between the
88078bcb3b5SPeter Brune    solution and the function.  For instance, SNESQN may be scaled by the
88178bcb3b5SPeter Brune    line search lambda using the argument -snes_qn_scaling ls.
88278bcb3b5SPeter Brune 
883f40b411bSPeter Brune 
884cd7522eaSPeter Brune .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetDamping(), SNESLineSearchApply()
885f40b411bSPeter Brune @*/
886f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda)
8876a388c36SPeter Brune {
8886a388c36SPeter Brune   PetscFunctionBegin;
889f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
8906a388c36SPeter Brune   PetscValidPointer(lambda, 2);
8916a388c36SPeter Brune   *lambda = linesearch->lambda;
8926a388c36SPeter Brune   PetscFunctionReturn(0);
8936a388c36SPeter Brune }
8946a388c36SPeter Brune 
8956a388c36SPeter Brune #undef __FUNCT__
896f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetLambda"
897f40b411bSPeter Brune /*@
898f1c6b773SPeter Brune    SNESLineSearchSetLambda - Sets the linesearch steplength.
899f40b411bSPeter Brune 
900f40b411bSPeter Brune    Input Parameters:
9018e557f58SPeter Brune +  linesearch - linesearch context
902f40b411bSPeter Brune -  lambda - The last steplength.
903f40b411bSPeter Brune 
904cd7522eaSPeter Brune    Notes:
905cd7522eaSPeter Brune    This routine is typically used within implementations of SNESLineSearchApply
906cd7522eaSPeter Brune    to set the final steplength.  This routine (and SNESLineSearchGetLambda()) were
907cd7522eaSPeter Brune    added in order to facilitate Quasi-Newton methods that use the previous steplength
908cd7522eaSPeter Brune    as an inner scaling parameter.
909cd7522eaSPeter Brune 
91078bcb3b5SPeter Brune    Level: advanced
911f40b411bSPeter Brune 
912f1c6b773SPeter Brune .seealso: SNESLineSearchGetLambda()
913f40b411bSPeter Brune @*/
914f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
9156a388c36SPeter Brune {
9166a388c36SPeter Brune   PetscFunctionBegin;
917f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
9186a388c36SPeter Brune   linesearch->lambda = lambda;
9196a388c36SPeter Brune   PetscFunctionReturn(0);
9206a388c36SPeter Brune }
9216a388c36SPeter Brune 
9226a388c36SPeter Brune #undef  __FUNCT__
923f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetTolerances"
924f40b411bSPeter Brune /*@
925*3c7d6663SPeter Brune    SNESLineSearchGetTolerances - Gets the tolerances for the linesearch.  These include
92678bcb3b5SPeter Brune    tolerances for the relative and absolute change in the function norm, the change
92778bcb3b5SPeter Brune    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
92878bcb3b5SPeter Brune    and the maximum number of iterations the line search procedure may take.
929f40b411bSPeter Brune 
930f40b411bSPeter Brune    Input Parameters:
9318e557f58SPeter Brune .  linesearch - linesearch context
932f40b411bSPeter Brune 
933f40b411bSPeter Brune    Output Parameters:
934516fe3c3SPeter Brune +  steptol - The minimum steplength
9356cc8e53bSPeter Brune .  maxstep - The maximum steplength
936516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
937516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
938516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
939516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
940f40b411bSPeter Brune 
94178bcb3b5SPeter Brune    Level: intermediate
94278bcb3b5SPeter Brune 
94378bcb3b5SPeter Brune    Notes:
94478bcb3b5SPeter Brune    Different line searches may implement these parameters slightly differently as
945*3c7d6663SPeter Brune    the type requires.
946516fe3c3SPeter Brune 
947f1c6b773SPeter Brune .seealso: SNESLineSearchSetTolerances()
948f40b411bSPeter Brune @*/
949f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
9506a388c36SPeter Brune {
9516a388c36SPeter Brune   PetscFunctionBegin;
952f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
953516fe3c3SPeter Brune   if (steptol) {
9546a388c36SPeter Brune     PetscValidPointer(steptol, 2);
9556a388c36SPeter Brune     *steptol = linesearch->steptol;
956516fe3c3SPeter Brune   }
957516fe3c3SPeter Brune   if (maxstep) {
958516fe3c3SPeter Brune     PetscValidPointer(maxstep, 3);
959516fe3c3SPeter Brune     *maxstep = linesearch->maxstep;
960516fe3c3SPeter Brune   }
961516fe3c3SPeter Brune   if (rtol) {
962516fe3c3SPeter Brune     PetscValidPointer(rtol, 4);
963516fe3c3SPeter Brune     *rtol = linesearch->rtol;
964516fe3c3SPeter Brune   }
965516fe3c3SPeter Brune   if (atol) {
966516fe3c3SPeter Brune     PetscValidPointer(atol, 5);
967516fe3c3SPeter Brune     *atol = linesearch->atol;
968516fe3c3SPeter Brune   }
969516fe3c3SPeter Brune   if (ltol) {
970516fe3c3SPeter Brune     PetscValidPointer(ltol, 6);
971516fe3c3SPeter Brune     *ltol = linesearch->ltol;
972516fe3c3SPeter Brune   }
973516fe3c3SPeter Brune   if (max_its) {
974516fe3c3SPeter Brune     PetscValidPointer(max_its, 7);
975516fe3c3SPeter Brune     *max_its = linesearch->max_its;
976516fe3c3SPeter Brune   }
9776a388c36SPeter Brune   PetscFunctionReturn(0);
9786a388c36SPeter Brune }
9796a388c36SPeter Brune 
9806a388c36SPeter Brune #undef  __FUNCT__
981f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetTolerances"
982f40b411bSPeter Brune /*@
983*3c7d6663SPeter Brune    SNESLineSearchSetTolerances -  Gets the tolerances for the linesearch.  These include
98478bcb3b5SPeter Brune    tolerances for the relative and absolute change in the function norm, the change
98578bcb3b5SPeter Brune    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
98678bcb3b5SPeter Brune    and the maximum number of iterations the line search procedure may take.
987f40b411bSPeter Brune 
988f40b411bSPeter Brune    Input Parameters:
9898e557f58SPeter Brune +  linesearch - linesearch context
990516fe3c3SPeter Brune .  steptol - The minimum steplength
9916cc8e53bSPeter Brune .  maxstep - The maximum steplength
992516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
993516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
994516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
995516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
996f40b411bSPeter Brune 
99778bcb3b5SPeter Brune    Notes:
998*3c7d6663SPeter Brune    The user may choose to not set any of the tolerances using PETSC_DEFAULT in
99978bcb3b5SPeter Brune    place of an argument.
1000f40b411bSPeter Brune 
100178bcb3b5SPeter Brune    Level: intermediate
1002516fe3c3SPeter Brune 
1003f1c6b773SPeter Brune .seealso: SNESLineSearchGetTolerances()
1004f40b411bSPeter Brune @*/
1005f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its)
10066a388c36SPeter Brune {
10076a388c36SPeter Brune   PetscFunctionBegin;
1008f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1009d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,steptol,2);
1010d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,maxstep,3);
1011d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,rtol,4);
1012d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,atol,5);
1013d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,ltol,6);
1014d3952184SSatish Balay   PetscValidLogicalCollectiveInt(linesearch,max_its,7);
1015d3952184SSatish Balay 
1016d3952184SSatish Balay   if ( steptol!= PETSC_DEFAULT) {
1017d3952184SSatish Balay     if ( steptol < 0.0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Minimum step length %G must be non-negative",steptol);
10186a388c36SPeter Brune     linesearch->steptol = steptol;
1019d3952184SSatish Balay   }
1020d3952184SSatish Balay 
1021d3952184SSatish Balay   if ( maxstep!= PETSC_DEFAULT) {
1022d3952184SSatish Balay     if ( maxstep < 0.0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Maximum step length %G must be non-negative",maxstep);
1023516fe3c3SPeter Brune     linesearch->maxstep = maxstep;
1024d3952184SSatish Balay   }
1025d3952184SSatish Balay 
1026d3952184SSatish Balay   if (rtol != PETSC_DEFAULT) {
1027d3952184SSatish 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);
1028516fe3c3SPeter Brune     linesearch->rtol = rtol;
1029d3952184SSatish Balay   }
1030d3952184SSatish Balay 
1031d3952184SSatish Balay   if (atol != PETSC_DEFAULT) {
1032d3952184SSatish Balay     if (atol < 0.0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %G must be non-negative",atol);
1033516fe3c3SPeter Brune     linesearch->atol = atol;
1034d3952184SSatish Balay   }
1035d3952184SSatish Balay 
1036d3952184SSatish Balay   if (ltol != PETSC_DEFAULT) {
1037d3952184SSatish Balay     if (ltol < 0.0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Labmda tolerance %G must be non-negative",ltol);
1038516fe3c3SPeter Brune   linesearch->ltol = ltol;
1039d3952184SSatish Balay   }
1040d3952184SSatish Balay 
1041d3952184SSatish Balay   if (max_its != PETSC_DEFAULT) {
1042d3952184SSatish Balay     if (max_its < 0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",max_its);
1043516fe3c3SPeter Brune     linesearch->max_its = max_its;
1044d3952184SSatish Balay   }
1045d3952184SSatish Balay 
10466a388c36SPeter Brune   PetscFunctionReturn(0);
10476a388c36SPeter Brune }
10486a388c36SPeter Brune 
1049516fe3c3SPeter Brune 
10506a388c36SPeter Brune #undef __FUNCT__
1051f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetDamping"
1052f40b411bSPeter Brune /*@
1053f1c6b773SPeter Brune    SNESLineSearchGetDamping - Gets the line search damping parameter.
1054f40b411bSPeter Brune 
1055f40b411bSPeter Brune    Input Parameters:
10568e557f58SPeter Brune .  linesearch - linesearch context
1057f40b411bSPeter Brune 
1058f40b411bSPeter Brune    Output Parameters:
10598e557f58SPeter Brune .  damping - The damping parameter
1060f40b411bSPeter Brune 
106178bcb3b5SPeter Brune    Level: advanced
1062f40b411bSPeter Brune 
106378bcb3b5SPeter Brune .seealso: SNESLineSearchGetStepTolerance(), SNESQN
1064f40b411bSPeter Brune @*/
1065f40b411bSPeter Brune 
1066f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping)
10676a388c36SPeter Brune {
10686a388c36SPeter Brune   PetscFunctionBegin;
1069f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
10706a388c36SPeter Brune   PetscValidPointer(damping, 2);
10716a388c36SPeter Brune   *damping = linesearch->damping;
10726a388c36SPeter Brune   PetscFunctionReturn(0);
10736a388c36SPeter Brune }
10746a388c36SPeter Brune 
10756a388c36SPeter Brune #undef __FUNCT__
1076f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetDamping"
1077f40b411bSPeter Brune /*@
1078f1c6b773SPeter Brune    SNESLineSearchSetDamping - Sets the line search damping paramter.
1079f40b411bSPeter Brune 
1080f40b411bSPeter Brune    Input Parameters:
108178bcb3b5SPeter Brune .  linesearch - linesearch context
108278bcb3b5SPeter Brune .  damping - The damping parameter
1083f40b411bSPeter Brune 
1084f40b411bSPeter Brune    Level: intermediate
1085f40b411bSPeter Brune 
1086cd7522eaSPeter Brune    Notes:
1087cd7522eaSPeter Brune    The basic line search merely takes the update step scaled by the damping parameter.
1088cd7522eaSPeter Brune    The use of the damping parameter in the l2 and cp line searches is much more subtle;
108978bcb3b5SPeter Brune    it is used as a starting point in calculating the secant step. However, the eventual
1090cd7522eaSPeter Brune    step may be of greater length than the damping parameter.  In the bt line search it is
1091cd7522eaSPeter Brune    used as the maximum possible step length, as the bt line search only backtracks.
1092cd7522eaSPeter Brune 
1093f1c6b773SPeter Brune .seealso: SNESLineSearchGetDamping()
1094f40b411bSPeter Brune @*/
1095f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping)
10966a388c36SPeter Brune {
10976a388c36SPeter Brune   PetscFunctionBegin;
1098f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
10996a388c36SPeter Brune   linesearch->damping = damping;
11006a388c36SPeter Brune   PetscFunctionReturn(0);
11016a388c36SPeter Brune }
11026a388c36SPeter Brune 
11036a388c36SPeter Brune #undef __FUNCT__
110459405d5eSPeter Brune #define __FUNCT__ "SNESLineSearchGetOrder"
110559405d5eSPeter Brune /*@
110659405d5eSPeter Brune    SNESLineSearchGetOrder - Gets the line search approximation order.
110759405d5eSPeter Brune 
110859405d5eSPeter Brune    Input Parameters:
110978bcb3b5SPeter Brune .  linesearch - linesearch context
111059405d5eSPeter Brune 
111159405d5eSPeter Brune    Output Parameters:
11128e557f58SPeter Brune .  order - The order
111359405d5eSPeter Brune 
111478bcb3b5SPeter Brune    Possible Values for order:
1115*3c7d6663SPeter Brune +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
1116*3c7d6663SPeter Brune .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
1117*3c7d6663SPeter Brune -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
111878bcb3b5SPeter Brune 
111959405d5eSPeter Brune    Level: intermediate
112059405d5eSPeter Brune 
112159405d5eSPeter Brune .seealso: SNESLineSearchSetOrder()
112259405d5eSPeter Brune @*/
112359405d5eSPeter Brune 
1124b000cd8dSPeter Brune PetscErrorCode  SNESLineSearchGetOrder(SNESLineSearch linesearch,PetscInt *order)
112559405d5eSPeter Brune {
112659405d5eSPeter Brune   PetscFunctionBegin;
112759405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
112859405d5eSPeter Brune   PetscValidPointer(order, 2);
112959405d5eSPeter Brune   *order = linesearch->order;
113059405d5eSPeter Brune   PetscFunctionReturn(0);
113159405d5eSPeter Brune }
113259405d5eSPeter Brune 
113359405d5eSPeter Brune #undef __FUNCT__
113459405d5eSPeter Brune #define __FUNCT__ "SNESLineSearchSetOrder"
113559405d5eSPeter Brune /*@
113659405d5eSPeter Brune    SNESLineSearchSetOrder - Sets the line search damping paramter.
113759405d5eSPeter Brune 
113859405d5eSPeter Brune    Input Parameters:
113978bcb3b5SPeter Brune .  linesearch - linesearch context
114078bcb3b5SPeter Brune .  order - The damping parameter
114159405d5eSPeter Brune 
114259405d5eSPeter Brune    Level: intermediate
114359405d5eSPeter Brune 
114478bcb3b5SPeter Brune    Possible Values for order:
1145*3c7d6663SPeter Brune +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
1146*3c7d6663SPeter Brune .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
1147*3c7d6663SPeter Brune -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
114878bcb3b5SPeter Brune 
114959405d5eSPeter Brune    Notes:
115059405d5eSPeter Brune    Variable orders are supported by the following line searches:
115178bcb3b5SPeter Brune +  bt - cubic and quadratic
115278bcb3b5SPeter Brune -  cp - linear and quadratic
115359405d5eSPeter Brune 
115459405d5eSPeter Brune .seealso: SNESLineSearchGetOrder()
115559405d5eSPeter Brune @*/
1156b000cd8dSPeter Brune PetscErrorCode  SNESLineSearchSetOrder(SNESLineSearch linesearch,PetscInt order)
115759405d5eSPeter Brune {
115859405d5eSPeter Brune   PetscFunctionBegin;
115959405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
116059405d5eSPeter Brune   linesearch->order = order;
116159405d5eSPeter Brune   PetscFunctionReturn(0);
116259405d5eSPeter Brune }
116359405d5eSPeter Brune 
116459405d5eSPeter Brune #undef __FUNCT__
1165f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetNorms"
1166f40b411bSPeter Brune /*@
1167f1c6b773SPeter Brune    SNESLineSearchGetNorms - Gets the norms for for X, Y, and F.
1168f40b411bSPeter Brune 
1169f40b411bSPeter Brune    Input Parameters:
117078bcb3b5SPeter Brune .  linesearch - linesearch context
1171f40b411bSPeter Brune 
1172f40b411bSPeter Brune    Output Parameters:
1173f40b411bSPeter Brune +  xnorm - The norm of the current solution
1174f40b411bSPeter Brune .  fnorm - The norm of the current function
1175f40b411bSPeter Brune -  ynorm - The norm of the current update
1176f40b411bSPeter Brune 
1177cd7522eaSPeter Brune    Notes:
1178cd7522eaSPeter Brune    This function is mainly called from SNES implementations.
1179cd7522eaSPeter Brune 
118078bcb3b5SPeter Brune    Level: developer
1181f40b411bSPeter Brune 
1182f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms() SNESLineSearchGetVecs()
1183f40b411bSPeter Brune @*/
1184f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm)
1185bf7f4e0aSPeter Brune {
1186bf7f4e0aSPeter Brune   PetscFunctionBegin;
1187f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1188bf7f4e0aSPeter Brune   if (xnorm) {
1189bf7f4e0aSPeter Brune     *xnorm = linesearch->xnorm;
1190bf7f4e0aSPeter Brune   }
1191bf7f4e0aSPeter Brune   if (fnorm) {
1192bf7f4e0aSPeter Brune     *fnorm = linesearch->fnorm;
1193bf7f4e0aSPeter Brune   }
1194bf7f4e0aSPeter Brune   if (ynorm) {
1195bf7f4e0aSPeter Brune     *ynorm = linesearch->ynorm;
1196bf7f4e0aSPeter Brune   }
1197bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1198bf7f4e0aSPeter Brune }
1199bf7f4e0aSPeter Brune 
1200e7058c64SPeter Brune #undef __FUNCT__
1201f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetNorms"
1202f40b411bSPeter Brune /*@
1203f1c6b773SPeter Brune    SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F.
1204f40b411bSPeter Brune 
1205f40b411bSPeter Brune    Input Parameters:
120678bcb3b5SPeter Brune +  linesearch - linesearch context
1207f40b411bSPeter Brune .  xnorm - The norm of the current solution
1208f40b411bSPeter Brune .  fnorm - The norm of the current function
1209f40b411bSPeter Brune -  ynorm - The norm of the current update
1210f40b411bSPeter Brune 
121178bcb3b5SPeter Brune    Level: advanced
1212f40b411bSPeter Brune 
1213f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1214f40b411bSPeter Brune @*/
1215f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
12166a388c36SPeter Brune {
12176a388c36SPeter Brune   PetscFunctionBegin;
1218f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
12196a388c36SPeter Brune   linesearch->xnorm = xnorm;
12206a388c36SPeter Brune   linesearch->fnorm = fnorm;
12216a388c36SPeter Brune   linesearch->ynorm = ynorm;
12226a388c36SPeter Brune   PetscFunctionReturn(0);
12236a388c36SPeter Brune }
12246a388c36SPeter Brune 
12256a388c36SPeter Brune #undef __FUNCT__
1226f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchComputeNorms"
1227f40b411bSPeter Brune /*@
1228f1c6b773SPeter Brune    SNESLineSearchComputeNorms - Computes the norms of X, F, and Y.
1229f40b411bSPeter Brune 
1230f40b411bSPeter Brune    Input Parameters:
123178bcb3b5SPeter Brune .  linesearch - linesearch context
1232f40b411bSPeter Brune 
1233f40b411bSPeter Brune    Options Database Keys:
12348e557f58SPeter Brune .   -snes_linesearch_norms - turn norm computation on or off
1235f40b411bSPeter Brune 
1236f40b411bSPeter Brune    Level: intermediate
1237f40b411bSPeter Brune 
123878bcb3b5SPeter Brune .seealso: SNESLineSearchGetNorms, SNESLineSearchSetNorms(), SNESLineSearchSetComputeNorms()
1239f40b411bSPeter Brune @*/
1240f1c6b773SPeter Brune PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
12416a388c36SPeter Brune {
12426a388c36SPeter Brune   PetscErrorCode ierr;
12439bd66eb0SPeter Brune   SNES snes;
12446a388c36SPeter Brune   PetscFunctionBegin;
12456a388c36SPeter Brune   if (linesearch->norms) {
12469bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
1247f1c6b773SPeter Brune       ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
12489bd66eb0SPeter Brune       ierr = VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
12499bd66eb0SPeter Brune       ierr = VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
12509bd66eb0SPeter Brune       ierr = (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);CHKERRQ(ierr);
12519bd66eb0SPeter Brune     } else {
12526a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_func,   NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
12536a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_sol,    NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
12546a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
12556a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_func,     NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
12566a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_sol,      NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
12576a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_update,   NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
12586a388c36SPeter Brune     }
12599bd66eb0SPeter Brune   }
12606a388c36SPeter Brune   PetscFunctionReturn(0);
12616a388c36SPeter Brune }
12626a388c36SPeter Brune 
12636f263ca3SPeter Brune 
12646f263ca3SPeter Brune #undef __FUNCT__
12656f263ca3SPeter Brune #define __FUNCT__ "SNESLineSearchSetComputeNorms"
12666f263ca3SPeter Brune /*@
12676f263ca3SPeter Brune    SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search.
12686f263ca3SPeter Brune 
12696f263ca3SPeter Brune    Input Parameters:
127078bcb3b5SPeter Brune +  linesearch  - linesearch context
127178bcb3b5SPeter Brune -  flg  - indicates whether or not to compute norms
12726f263ca3SPeter Brune 
12736f263ca3SPeter Brune    Options Database Keys:
12748e557f58SPeter Brune .   -snes_linesearch_norms - turn norm computation on or off
12756f263ca3SPeter Brune 
12766f263ca3SPeter Brune    Notes:
12771a4f838cSPeter Brune    This is most relevant to the SNESLINESEARCHBASIC line search type.
12786f263ca3SPeter Brune 
12796f263ca3SPeter Brune    Level: intermediate
12806f263ca3SPeter Brune 
12811a4f838cSPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetNorms(), SNESLineSearchComputeNorms(), SNESLINESEARCHBASIC
12826f263ca3SPeter Brune @*/
12836f263ca3SPeter Brune PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg)
12846f263ca3SPeter Brune {
12856f263ca3SPeter Brune   PetscFunctionBegin;
12866f263ca3SPeter Brune   linesearch->norms = flg;
12876f263ca3SPeter Brune   PetscFunctionReturn(0);
12886f263ca3SPeter Brune }
12896f263ca3SPeter Brune 
12906a388c36SPeter Brune #undef __FUNCT__
1291f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetVecs"
1292f40b411bSPeter Brune /*@
1293f1c6b773SPeter Brune    SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context
1294f40b411bSPeter Brune 
1295f40b411bSPeter Brune    Input Parameters:
129678bcb3b5SPeter Brune .  linesearch - linesearch context
1297f40b411bSPeter Brune 
1298f40b411bSPeter Brune    Output Parameters:
1299f40b411bSPeter Brune +  X - The old solution
1300f40b411bSPeter Brune .  F - The old function
1301f40b411bSPeter Brune .  Y - The search direction
1302f40b411bSPeter Brune .  W - The new solution
1303f40b411bSPeter Brune -  G - The new function
1304f40b411bSPeter Brune 
130578bcb3b5SPeter Brune    Level: advanced
1306f40b411bSPeter Brune 
1307f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1308f40b411bSPeter Brune @*/
1309f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G) {
13106a388c36SPeter Brune   PetscFunctionBegin;
1311f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
13126a388c36SPeter Brune   if (X) {
13136a388c36SPeter Brune     PetscValidPointer(X, 2);
13146a388c36SPeter Brune     *X = linesearch->vec_sol;
13156a388c36SPeter Brune   }
13166a388c36SPeter Brune   if (F) {
13176a388c36SPeter Brune     PetscValidPointer(F, 3);
13186a388c36SPeter Brune     *F = linesearch->vec_func;
13196a388c36SPeter Brune   }
13206a388c36SPeter Brune   if (Y) {
13216a388c36SPeter Brune     PetscValidPointer(Y, 4);
13226a388c36SPeter Brune     *Y = linesearch->vec_update;
13236a388c36SPeter Brune   }
13246a388c36SPeter Brune   if (W) {
13256a388c36SPeter Brune     PetscValidPointer(W, 5);
13266a388c36SPeter Brune     *W = linesearch->vec_sol_new;
13276a388c36SPeter Brune   }
13286a388c36SPeter Brune   if (G) {
13296a388c36SPeter Brune     PetscValidPointer(G, 6);
13306a388c36SPeter Brune     *G = linesearch->vec_func_new;
13316a388c36SPeter Brune   }
13326a388c36SPeter Brune 
13336a388c36SPeter Brune   PetscFunctionReturn(0);
13346a388c36SPeter Brune }
13356a388c36SPeter Brune 
13366a388c36SPeter Brune #undef __FUNCT__
1337f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetVecs"
1338f40b411bSPeter Brune /*@
1339f1c6b773SPeter Brune    SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context
1340f40b411bSPeter Brune 
1341f40b411bSPeter Brune    Input Parameters:
134278bcb3b5SPeter Brune +  linesearch - linesearch context
1343f40b411bSPeter Brune .  X - The old solution
1344f40b411bSPeter Brune .  F - The old function
1345f40b411bSPeter Brune .  Y - The search direction
1346f40b411bSPeter Brune .  W - The new solution
1347f40b411bSPeter Brune -  G - The new function
1348f40b411bSPeter Brune 
134978bcb3b5SPeter Brune    Level: advanced
1350f40b411bSPeter Brune 
1351f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms(), SNESLineSearchGetVecs()
1352f40b411bSPeter Brune @*/
1353f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G) {
13546a388c36SPeter Brune   PetscFunctionBegin;
1355f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
13566a388c36SPeter Brune   if (X) {
13576a388c36SPeter Brune     PetscValidHeaderSpecific(X,VEC_CLASSID,2);
13586a388c36SPeter Brune     linesearch->vec_sol = X;
13596a388c36SPeter Brune   }
13606a388c36SPeter Brune   if (F) {
13616a388c36SPeter Brune     PetscValidHeaderSpecific(F,VEC_CLASSID,3);
13626a388c36SPeter Brune     linesearch->vec_func = F;
13636a388c36SPeter Brune   }
13646a388c36SPeter Brune   if (Y) {
13656a388c36SPeter Brune     PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
13666a388c36SPeter Brune     linesearch->vec_update = Y;
13676a388c36SPeter Brune   }
13686a388c36SPeter Brune   if (W) {
13696a388c36SPeter Brune     PetscValidHeaderSpecific(W,VEC_CLASSID,5);
13706a388c36SPeter Brune     linesearch->vec_sol_new = W;
13716a388c36SPeter Brune   }
13726a388c36SPeter Brune   if (G) {
13736a388c36SPeter Brune     PetscValidHeaderSpecific(G,VEC_CLASSID,6);
13746a388c36SPeter Brune     linesearch->vec_func_new = G;
13756a388c36SPeter Brune   }
13766a388c36SPeter Brune 
13776a388c36SPeter Brune   PetscFunctionReturn(0);
13786a388c36SPeter Brune }
13796a388c36SPeter Brune 
13806a388c36SPeter Brune #undef __FUNCT__
1381f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchAppendOptionsPrefix"
1382e7058c64SPeter Brune /*@C
1383f1c6b773SPeter Brune    SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1384e7058c64SPeter Brune    SNES options in the database.
1385e7058c64SPeter Brune 
1386cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
1387e7058c64SPeter Brune 
1388e7058c64SPeter Brune    Input Parameters:
1389e7058c64SPeter Brune +  snes - the SNES context
1390e7058c64SPeter Brune -  prefix - the prefix to prepend to all option names
1391e7058c64SPeter Brune 
1392e7058c64SPeter Brune    Notes:
1393e7058c64SPeter Brune    A hyphen (-) must NOT be given at the beginning of the prefix name.
1394e7058c64SPeter Brune    The first character of all runtime options is AUTOMATICALLY the hyphen.
1395e7058c64SPeter Brune 
1396e7058c64SPeter Brune    Level: advanced
1397e7058c64SPeter Brune 
1398f1c6b773SPeter Brune .keywords: SNESLineSearch, append, options, prefix, database
1399e7058c64SPeter Brune 
1400e7058c64SPeter Brune .seealso: SNESGetOptionsPrefix()
1401e7058c64SPeter Brune @*/
1402f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[])
1403e7058c64SPeter Brune {
1404e7058c64SPeter Brune   PetscErrorCode ierr;
1405e7058c64SPeter Brune 
1406e7058c64SPeter Brune   PetscFunctionBegin;
1407f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1408e7058c64SPeter Brune   ierr = PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1409e7058c64SPeter Brune   PetscFunctionReturn(0);
1410e7058c64SPeter Brune }
1411e7058c64SPeter Brune 
1412e7058c64SPeter Brune #undef __FUNCT__
1413f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetOptionsPrefix"
1414e7058c64SPeter Brune /*@C
1415f1c6b773SPeter Brune    SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all
1416f1c6b773SPeter Brune    SNESLineSearch options in the database.
1417e7058c64SPeter Brune 
1418e7058c64SPeter Brune    Not Collective
1419e7058c64SPeter Brune 
1420e7058c64SPeter Brune    Input Parameter:
1421cd7522eaSPeter Brune .  linesearch - the SNESLineSearch context
1422e7058c64SPeter Brune 
1423e7058c64SPeter Brune    Output Parameter:
1424e7058c64SPeter Brune .  prefix - pointer to the prefix string used
1425e7058c64SPeter Brune 
14268e557f58SPeter Brune    Notes:
14278e557f58SPeter Brune    On the fortran side, the user should pass in a string 'prefix' of
1428e7058c64SPeter Brune    sufficient length to hold the prefix.
1429e7058c64SPeter Brune 
1430e7058c64SPeter Brune    Level: advanced
1431e7058c64SPeter Brune 
1432f1c6b773SPeter Brune .keywords: SNESLineSearch, get, options, prefix, database
1433e7058c64SPeter Brune 
1434e7058c64SPeter Brune .seealso: SNESAppendOptionsPrefix()
1435e7058c64SPeter Brune @*/
1436f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[])
1437e7058c64SPeter Brune {
1438e7058c64SPeter Brune   PetscErrorCode ierr;
1439e7058c64SPeter Brune 
1440e7058c64SPeter Brune   PetscFunctionBegin;
1441f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1442e7058c64SPeter Brune   ierr = PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1443e7058c64SPeter Brune   PetscFunctionReturn(0);
1444e7058c64SPeter Brune }
1445bf7f4e0aSPeter Brune 
1446bf7f4e0aSPeter Brune #undef __FUNCT__
1447f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetWork"
1448f40b411bSPeter Brune /*@
1449f1c6b773SPeter Brune    SNESLineSearchGetWork - Gets work vectors for the line search.
1450f40b411bSPeter Brune 
1451f40b411bSPeter Brune    Input Parameter:
1452f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
1453f40b411bSPeter Brune -  nwork - the number of work vectors
1454f40b411bSPeter Brune 
1455f40b411bSPeter Brune    Level: developer
1456f40b411bSPeter Brune 
1457cd7522eaSPeter Brune    Notes:
1458cd7522eaSPeter Brune    This is typically called at the beginning of a SNESLineSearch or SNESLineSearchShell implementation.
1459cd7522eaSPeter Brune 
1460f1c6b773SPeter Brune .keywords: SNESLineSearch, work, vector
1461f40b411bSPeter Brune 
1462f40b411bSPeter Brune .seealso: SNESDefaultGetWork()
1463f40b411bSPeter Brune @*/
1464f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetWork(SNESLineSearch linesearch, PetscInt nwork)
1465bf7f4e0aSPeter Brune {
1466bf7f4e0aSPeter Brune   PetscErrorCode ierr;
1467bf7f4e0aSPeter Brune   PetscFunctionBegin;
1468bf7f4e0aSPeter Brune   if (linesearch->vec_sol) {
1469bf7f4e0aSPeter Brune     ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr);
1470bf7f4e0aSPeter Brune   } else {
1471bf7f4e0aSPeter Brune     SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
1472bf7f4e0aSPeter Brune   }
1473bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1474bf7f4e0aSPeter Brune }
1475bf7f4e0aSPeter Brune 
14766a388c36SPeter Brune 
1477bf7f4e0aSPeter Brune #undef __FUNCT__
1478f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetSuccess"
1479f40b411bSPeter Brune /*@
1480f1c6b773SPeter Brune    SNESLineSearchGetSuccess - Gets the success/failure status of the last line search application
1481f40b411bSPeter Brune 
1482f40b411bSPeter Brune    Input Parameters:
148378bcb3b5SPeter Brune .  linesearch - linesearch context
1484f40b411bSPeter Brune 
1485f40b411bSPeter Brune    Output Parameters:
14868e557f58SPeter Brune .  success - The success or failure status
1487f40b411bSPeter Brune 
1488cd7522eaSPeter Brune    Notes:
1489cd7522eaSPeter Brune    This is typically called after SNESLineSearchApply in order to determine if the line-search failed
1490cd7522eaSPeter Brune    (and set the SNES convergence accordingly).
1491cd7522eaSPeter Brune 
1492f40b411bSPeter Brune    Level: intermediate
1493f40b411bSPeter Brune 
1494f1c6b773SPeter Brune .seealso: SNESLineSearchSetSuccess()
1495f40b411bSPeter Brune @*/
1496f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetSuccess(SNESLineSearch linesearch, PetscBool *success)
1497bf7f4e0aSPeter Brune {
1498bf7f4e0aSPeter Brune   PetscFunctionBegin;
1499f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
15006a388c36SPeter Brune   PetscValidPointer(success, 2);
1501bf7f4e0aSPeter Brune   if (success) {
1502bf7f4e0aSPeter Brune     *success = linesearch->success;
1503bf7f4e0aSPeter Brune   }
1504bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1505bf7f4e0aSPeter Brune }
1506bf7f4e0aSPeter Brune 
1507bf7f4e0aSPeter Brune #undef __FUNCT__
1508f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetSuccess"
1509f40b411bSPeter Brune /*@
1510f1c6b773SPeter Brune    SNESLineSearchSetSuccess - Sets the success/failure status of the last line search application
1511f40b411bSPeter Brune 
1512f40b411bSPeter Brune    Input Parameters:
151378bcb3b5SPeter Brune +  linesearch - linesearch context
15148e557f58SPeter Brune -  success - The success or failure status
1515f40b411bSPeter Brune 
1516cd7522eaSPeter Brune    Notes:
1517cd7522eaSPeter Brune    This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set
1518cd7522eaSPeter Brune    the success or failure of the line search method.
1519cd7522eaSPeter Brune 
152078bcb3b5SPeter Brune    Level: developer
1521f40b411bSPeter Brune 
1522f1c6b773SPeter Brune .seealso: SNESLineSearchGetSuccess()
1523f40b411bSPeter Brune @*/
1524f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetSuccess(SNESLineSearch linesearch, PetscBool success)
15256a388c36SPeter Brune {
1526f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
15276a388c36SPeter Brune   PetscFunctionBegin;
15286a388c36SPeter Brune   linesearch->success = success;
15296a388c36SPeter Brune   PetscFunctionReturn(0);
15306a388c36SPeter Brune }
15316a388c36SPeter Brune 
15326a388c36SPeter Brune #undef __FUNCT__
1533f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetVIFunctions"
15349bd66eb0SPeter Brune /*@C
1535f1c6b773SPeter Brune    SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.
15369bd66eb0SPeter Brune 
15379bd66eb0SPeter Brune    Input Parameters:
15389bd66eb0SPeter Brune +  snes - nonlinear context obtained from SNESCreate()
15399bd66eb0SPeter Brune .  projectfunc - function for projecting the function to the bounds
15409bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
15419bd66eb0SPeter Brune 
15429bd66eb0SPeter Brune    Logically Collective on SNES
15439bd66eb0SPeter Brune 
15449bd66eb0SPeter Brune    Calling sequence of projectfunc:
15459bd66eb0SPeter Brune .vb
15469bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X)
15479bd66eb0SPeter Brune .ve
15489bd66eb0SPeter Brune 
15499bd66eb0SPeter Brune     Input parameters for projectfunc:
15509bd66eb0SPeter Brune +   snes - nonlinear context
15519bd66eb0SPeter Brune -   X - current solution
15529bd66eb0SPeter Brune 
1553cd7522eaSPeter Brune     Output parameters for projectfunc:
15549bd66eb0SPeter Brune .   X - Projected solution
15559bd66eb0SPeter Brune 
15569bd66eb0SPeter Brune    Calling sequence of normfunc:
15579bd66eb0SPeter Brune .vb
15589bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm)
15599bd66eb0SPeter Brune .ve
15609bd66eb0SPeter Brune 
1561cd7522eaSPeter Brune     Input parameters for normfunc:
15629bd66eb0SPeter Brune +   snes - nonlinear context
15639bd66eb0SPeter Brune .   X - current solution
15649bd66eb0SPeter Brune -   F - current residual
15659bd66eb0SPeter Brune 
1566cd7522eaSPeter Brune     Output parameters for normfunc:
15679bd66eb0SPeter Brune .   fnorm - VI-specific norm of the function
15689bd66eb0SPeter Brune 
1569cd7522eaSPeter Brune     Notes:
1570cd7522eaSPeter Brune     The VI solvers require projection of the solution to the feasible set.  projectfunc should implement this.
1571cd7522eaSPeter Brune 
1572cd7522eaSPeter Brune     The VI solvers require special evaluation of the function norm such that the norm is only calculated
1573cd7522eaSPeter Brune     on the inactive set.  This should be implemented by normfunc.
15749bd66eb0SPeter Brune 
15759bd66eb0SPeter Brune     Level: developer
15769bd66eb0SPeter Brune 
15779bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, set, line search
15789bd66eb0SPeter Brune 
1579f1c6b773SPeter Brune .seealso: SNESLineSearchGetVIFunctions(), SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck()
15809bd66eb0SPeter Brune @*/
1581f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc)
15829bd66eb0SPeter Brune {
15839bd66eb0SPeter Brune   PetscFunctionBegin;
1584f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
15859bd66eb0SPeter Brune   if (projectfunc) linesearch->ops->viproject = projectfunc;
15869bd66eb0SPeter Brune   if (normfunc) linesearch->ops->vinorm = normfunc;
15879bd66eb0SPeter Brune   PetscFunctionReturn(0);
15889bd66eb0SPeter Brune }
15899bd66eb0SPeter Brune 
15909bd66eb0SPeter Brune #undef __FUNCT__
1591f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetVIFunctions"
15929bd66eb0SPeter Brune /*@C
1593f1c6b773SPeter Brune    SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.
15949bd66eb0SPeter Brune 
15959bd66eb0SPeter Brune    Input Parameters:
15969bd66eb0SPeter Brune .  snes - nonlinear context obtained from SNESCreate()
15979bd66eb0SPeter Brune 
15989bd66eb0SPeter Brune    Output Parameters:
15999bd66eb0SPeter Brune +  projectfunc - function for projecting the function to the bounds
16009bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
16019bd66eb0SPeter Brune 
16029bd66eb0SPeter Brune    Logically Collective on SNES
16039bd66eb0SPeter Brune 
16049bd66eb0SPeter Brune     Level: developer
16059bd66eb0SPeter Brune 
16069bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, get, line search
16079bd66eb0SPeter Brune 
1608f1c6b773SPeter Brune .seealso: SNESLineSearchSetVIFunctions(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
16099bd66eb0SPeter Brune @*/
1610f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc)
16119bd66eb0SPeter Brune {
16129bd66eb0SPeter Brune   PetscFunctionBegin;
16139bd66eb0SPeter Brune   if (projectfunc) *projectfunc = linesearch->ops->viproject;
16149bd66eb0SPeter Brune   if (normfunc) *normfunc = linesearch->ops->vinorm;
16159bd66eb0SPeter Brune   PetscFunctionReturn(0);
16169bd66eb0SPeter Brune }
16179bd66eb0SPeter Brune 
16189bd66eb0SPeter Brune #undef __FUNCT__
1619f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchRegister"
1620bf7f4e0aSPeter Brune /*@C
1621f1c6b773SPeter Brune   SNESLineSearchRegister - See SNESLineSearchRegisterDynamic()
1622bf7f4e0aSPeter Brune 
1623bf7f4e0aSPeter Brune   Level: advanced
1624bf7f4e0aSPeter Brune @*/
1625f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNESLineSearch))
1626bf7f4e0aSPeter Brune {
1627bf7f4e0aSPeter Brune   char           fullname[PETSC_MAX_PATH_LEN];
1628bf7f4e0aSPeter Brune   PetscErrorCode ierr;
1629bf7f4e0aSPeter Brune 
1630bf7f4e0aSPeter Brune   PetscFunctionBegin;
1631bf7f4e0aSPeter Brune   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
1632f1c6b773SPeter Brune   ierr = PetscFListAdd(&SNESLineSearchList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
1633bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1634bf7f4e0aSPeter Brune }
1635