xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision 38bcdd5a967c22b143d2cc59677183600d4e9763)
1b45d2f2cSJed Brown #include <petsc-private/linesearchimpl.h> /*I "petscsnes.h" I*/
2bf7f4e0aSPeter Brune 
3f1c6b773SPeter Brune PetscBool  SNESLineSearchRegisterAllCalled = PETSC_FALSE;
4f1c6b773SPeter Brune PetscFList SNESLineSearchList              = PETSC_NULL;
5bf7f4e0aSPeter Brune 
6f1c6b773SPeter Brune PetscClassId   SNESLINESEARCH_CLASSID;
7f1c6b773SPeter Brune PetscLogEvent  SNESLineSearch_Apply;
8bf7f4e0aSPeter Brune 
9bf7f4e0aSPeter Brune #undef __FUNCT__
10f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate"
11f40b411bSPeter Brune /*@
12cd7522eaSPeter Brune    SNESLineSearchCreate - Creates the line search context.
13f40b411bSPeter Brune 
14cd7522eaSPeter Brune    Logically Collective on Comm
15f40b411bSPeter Brune 
16f40b411bSPeter Brune    Input Parameters:
17cd7522eaSPeter Brune .  comm - MPI communicator for the line search (typically from the associated SNES context).
18f40b411bSPeter Brune 
19f40b411bSPeter Brune    Output Parameters:
208e557f58SPeter Brune .  outlinesearch - the new linesearch context
21f40b411bSPeter Brune 
2278bcb3b5SPeter Brune    Level: beginner
23f40b411bSPeter Brune 
24cd7522eaSPeter Brune .keywords: LineSearch, create, context
25f40b411bSPeter Brune 
26f40b411bSPeter Brune .seealso: LineSearchDestroy()
27f40b411bSPeter Brune @*/
28f40b411bSPeter Brune 
29f1c6b773SPeter Brune PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch) {
30bf7f4e0aSPeter Brune   PetscErrorCode      ierr;
31f1c6b773SPeter Brune   SNESLineSearch     linesearch;
32bf7f4e0aSPeter Brune   PetscFunctionBegin;
33ea5d4fccSPeter Brune   PetscValidPointer(outlinesearch,2);
34ea5d4fccSPeter Brune   *outlinesearch = PETSC_NULL;
35f1c6b773SPeter Brune   ierr = PetscHeaderCreate(linesearch,_p_LineSearch,struct _LineSearchOps,SNESLINESEARCH_CLASSID, 0,
36f1c6b773SPeter Brune                            "SNESLineSearch","Line-search method","SNESLineSearch",comm,SNESLineSearchDestroy,SNESLineSearchView);CHKERRQ(ierr);
37bf7f4e0aSPeter Brune 
38bf7f4e0aSPeter Brune   linesearch->ops->precheckstep = PETSC_NULL;
39bf7f4e0aSPeter Brune   linesearch->ops->postcheckstep = PETSC_NULL;
40bf7f4e0aSPeter Brune 
419bd66eb0SPeter Brune   linesearch->vec_sol_new   = PETSC_NULL;
429bd66eb0SPeter Brune   linesearch->vec_func_new  = PETSC_NULL;
439bd66eb0SPeter Brune   linesearch->vec_sol       = PETSC_NULL;
449bd66eb0SPeter Brune   linesearch->vec_func      = PETSC_NULL;
459bd66eb0SPeter Brune   linesearch->vec_update    = PETSC_NULL;
469bd66eb0SPeter Brune 
47bf7f4e0aSPeter Brune   linesearch->lambda        = 1.0;
48bf7f4e0aSPeter Brune   linesearch->fnorm         = 1.0;
49bf7f4e0aSPeter Brune   linesearch->ynorm         = 1.0;
50bf7f4e0aSPeter Brune   linesearch->xnorm         = 1.0;
51bf7f4e0aSPeter Brune   linesearch->success       = PETSC_TRUE;
52bf7f4e0aSPeter Brune   linesearch->norms         = PETSC_TRUE;
53bf7f4e0aSPeter Brune   linesearch->keeplambda    = PETSC_FALSE;
54bf7f4e0aSPeter Brune   linesearch->damping       = 1.0;
55bf7f4e0aSPeter Brune   linesearch->maxstep       = 1e8;
56bf7f4e0aSPeter Brune   linesearch->steptol       = 1e-12;
57516fe3c3SPeter Brune   linesearch->rtol          = 1e-8;
58516fe3c3SPeter Brune   linesearch->atol          = 1e-15;
59516fe3c3SPeter Brune   linesearch->ltol          = 1e-8;
60bf7f4e0aSPeter Brune   linesearch->precheckctx   = PETSC_NULL;
61bf7f4e0aSPeter Brune   linesearch->postcheckctx  = PETSC_NULL;
6259405d5eSPeter Brune   linesearch->max_its       = 1;
63bf7f4e0aSPeter Brune   linesearch->setupcalled   = PETSC_FALSE;
64bf7f4e0aSPeter Brune   *outlinesearch            = linesearch;
65bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
66bf7f4e0aSPeter Brune }
67bf7f4e0aSPeter Brune 
68bf7f4e0aSPeter Brune #undef __FUNCT__
69f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetUp"
70f40b411bSPeter Brune /*@
7178bcb3b5SPeter Brune    SNESLineSearchSetUp - Prepares the line search for being applied by allocating
7278bcb3b5SPeter Brune    any required vectors.
73f40b411bSPeter Brune 
74cd7522eaSPeter Brune    Collective on SNESLineSearch
75f40b411bSPeter Brune 
76f40b411bSPeter Brune    Input Parameters:
77f40b411bSPeter Brune .  linesearch - The LineSearch instance.
78f40b411bSPeter Brune 
79cd7522eaSPeter Brune    Notes:
80cd7522eaSPeter Brune    For most cases, this needn't be called outside of SNESLineSearchApply().
81cd7522eaSPeter Brune    The only current case where this is called outside of this is for the VI
8278bcb3b5SPeter Brune    solvers, which modify the solution and work vectors before the first call
83cd7522eaSPeter Brune    of SNESLineSearchApply, requiring the SNESLineSearch work vectors to be
84cd7522eaSPeter Brune    allocated upfront.
85cd7522eaSPeter Brune 
86cd7522eaSPeter Brune 
8778bcb3b5SPeter Brune    Level: advanced
88f40b411bSPeter Brune 
89f1c6b773SPeter Brune .keywords: SNESLineSearch, SetUp
90f40b411bSPeter Brune 
91f1c6b773SPeter Brune .seealso: SNESLineSearchReset()
92f40b411bSPeter Brune @*/
93f40b411bSPeter Brune 
94f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch) {
95bf7f4e0aSPeter Brune   PetscErrorCode ierr;
96bf7f4e0aSPeter Brune   PetscFunctionBegin;
97bf7f4e0aSPeter Brune   if (!((PetscObject)linesearch)->type_name) {
981a4f838cSPeter Brune     ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr);
99bf7f4e0aSPeter Brune   }
100bf7f4e0aSPeter Brune   if (!linesearch->setupcalled) {
1019bd66eb0SPeter Brune     if (!linesearch->vec_sol_new) {
102bf7f4e0aSPeter Brune       ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new);CHKERRQ(ierr);
1039bd66eb0SPeter Brune     }
1049bd66eb0SPeter Brune     if (!linesearch->vec_func_new) {
1059bd66eb0SPeter Brune       ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new);CHKERRQ(ierr);
1069bd66eb0SPeter Brune     }
107bf7f4e0aSPeter Brune     if (linesearch->ops->setup) {
108bf7f4e0aSPeter Brune       ierr = (*linesearch->ops->setup)(linesearch);CHKERRQ(ierr);
109bf7f4e0aSPeter Brune     }
110bf7f4e0aSPeter Brune     linesearch->lambda = linesearch->damping;
111bf7f4e0aSPeter Brune     linesearch->setupcalled = PETSC_TRUE;
112bf7f4e0aSPeter Brune   }
113bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
114bf7f4e0aSPeter Brune }
115bf7f4e0aSPeter Brune 
116bf7f4e0aSPeter Brune #undef __FUNCT__
117f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchReset"
118f40b411bSPeter Brune 
119f40b411bSPeter Brune /*@
12078bcb3b5SPeter Brune    SNESLineSearchReset - Undoes the SetUp and deletes any Vecs or Mats allocated by the line search.
121f40b411bSPeter Brune 
122f1c6b773SPeter Brune    Collective on SNESLineSearch
123f40b411bSPeter Brune 
124f40b411bSPeter Brune    Input Parameters:
125f40b411bSPeter Brune .  linesearch - The LineSearch instance.
126f40b411bSPeter Brune 
12778bcb3b5SPeter Brune    Level: intermediate
128f40b411bSPeter Brune 
129cd7522eaSPeter Brune .keywords: SNESLineSearch, Reset
130f40b411bSPeter Brune 
131f1c6b773SPeter Brune .seealso: SNESLineSearchSetUp()
132f40b411bSPeter Brune @*/
133f40b411bSPeter Brune 
134f1c6b773SPeter Brune PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch) {
135bf7f4e0aSPeter Brune   PetscErrorCode ierr;
136bf7f4e0aSPeter Brune   PetscFunctionBegin;
137bf7f4e0aSPeter Brune   if (linesearch->ops->reset) {
138bf7f4e0aSPeter Brune     (*linesearch->ops->reset)(linesearch);
139bf7f4e0aSPeter Brune   }
140bf7f4e0aSPeter Brune   ierr = VecDestroy(&linesearch->vec_sol_new);CHKERRQ(ierr);
141bf7f4e0aSPeter Brune   ierr = VecDestroy(&linesearch->vec_func_new);CHKERRQ(ierr);
142bf7f4e0aSPeter Brune 
143bf7f4e0aSPeter Brune   ierr = VecDestroyVecs(linesearch->nwork, &linesearch->work);CHKERRQ(ierr);
144bf7f4e0aSPeter Brune   linesearch->nwork = 0;
145bf7f4e0aSPeter Brune   linesearch->setupcalled = PETSC_FALSE;
146bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
147bf7f4e0aSPeter Brune }
148bf7f4e0aSPeter Brune 
14986d74e61SPeter Brune 
15086d74e61SPeter Brune #undef __FUNCT__
151f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetPreCheck"
15286d74e61SPeter Brune /*@C
153f1c6b773SPeter Brune    SNESLineSearchSetPreCheck - Sets a pre-check function for the line search routine.
15486d74e61SPeter Brune 
155f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
15686d74e61SPeter Brune 
15786d74e61SPeter Brune    Input Parameters:
158f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
15986d74e61SPeter Brune .  func       - [optional] function evaluation routine
16086d74e61SPeter Brune -  ctx        - [optional] user-defined context for private data for the
16186d74e61SPeter Brune                 function evaluation routine (may be PETSC_NULL)
16286d74e61SPeter Brune 
16386d74e61SPeter Brune    Calling sequence of func:
164f1c6b773SPeter Brune $    func (SNESLineSearch snes,Vec x,Vec y, PetscBool *changed);
16586d74e61SPeter Brune 
16686d74e61SPeter Brune +  x - solution vector
16786d74e61SPeter Brune .  y - search direction vector
168cd7522eaSPeter Brune -  changed - flag to indicate the precheck changed x or y.
16986d74e61SPeter Brune 
17086d74e61SPeter Brune    Level: intermediate
17186d74e61SPeter Brune 
17286d74e61SPeter Brune .keywords: set, linesearch, pre-check
17386d74e61SPeter Brune 
174f1c6b773SPeter Brune .seealso: SNESLineSearchSetPostCheck()
17586d74e61SPeter Brune @*/
176f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetPreCheck(SNESLineSearch linesearch, SNESLineSearchPreCheckFunc func,void *ctx)
17786d74e61SPeter Brune {
1789bd66eb0SPeter Brune   PetscFunctionBegin;
179f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
18086d74e61SPeter Brune   if (func) linesearch->ops->precheckstep = func;
18186d74e61SPeter Brune   if (ctx) linesearch->precheckctx = ctx;
18286d74e61SPeter Brune   PetscFunctionReturn(0);
18386d74e61SPeter Brune }
18486d74e61SPeter Brune 
18586d74e61SPeter Brune 
18686d74e61SPeter Brune #undef __FUNCT__
187f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetPreCheck"
18886d74e61SPeter Brune /*@C
189cd7522eaSPeter Brune    SNESLineSearchSetPreCheck - Gets the pre-check function for the line search routine.
19086d74e61SPeter Brune 
19186d74e61SPeter Brune    Input Parameters:
192f1c6b773SPeter Brune .  linesearch - the SNESLineSearch context
19386d74e61SPeter Brune 
19486d74e61SPeter Brune    Output Parameters:
19586d74e61SPeter Brune +  func       - [optional] function evaluation routine
19686d74e61SPeter Brune -  ctx        - [optional] user-defined context for private data for the
19786d74e61SPeter Brune                 function evaluation routine (may be PETSC_NULL)
19886d74e61SPeter Brune 
19986d74e61SPeter Brune    Level: intermediate
20086d74e61SPeter Brune 
20186d74e61SPeter Brune .keywords: get, linesearch, pre-check
20286d74e61SPeter Brune 
203f1c6b773SPeter Brune .seealso: SNESLineSearchGetPostCheck(), SNESLineSearchSetPreCheck()
20486d74e61SPeter Brune @*/
205f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetPreCheck(SNESLineSearch linesearch, SNESLineSearchPreCheckFunc *func,void **ctx)
20686d74e61SPeter Brune {
2079bd66eb0SPeter Brune   PetscFunctionBegin;
208f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
20986d74e61SPeter Brune   if (func) *func = linesearch->ops->precheckstep;
21086d74e61SPeter Brune   if (ctx) *ctx = linesearch->precheckctx;
21186d74e61SPeter Brune   PetscFunctionReturn(0);
21286d74e61SPeter Brune }
21386d74e61SPeter Brune 
21486d74e61SPeter Brune 
21586d74e61SPeter Brune #undef __FUNCT__
216f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetPostCheck"
21786d74e61SPeter Brune /*@C
218f1c6b773SPeter Brune    SNESLineSearchSetPostCheck - Sets a post-check function for the line search routine.
21986d74e61SPeter Brune 
220f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
22186d74e61SPeter Brune 
22286d74e61SPeter Brune    Input Parameters:
223f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
22486d74e61SPeter Brune .  func       - [optional] function evaluation routine
22586d74e61SPeter Brune -  ctx        - [optional] user-defined context for private data for the
22686d74e61SPeter Brune                 function evaluation routine (may be PETSC_NULL)
22786d74e61SPeter Brune 
22886d74e61SPeter Brune    Calling sequence of func:
229f1c6b773SPeter Brune $    func (SNESLineSearch linesearch,Vec x, Vec w, Vec y, PetscBool *changed_w, *changed_y);
23086d74e61SPeter Brune 
23186d74e61SPeter Brune +  x - old solution vector
23286d74e61SPeter Brune .  y - search direction vector
23386d74e61SPeter Brune .  w - new solution vector
2348e557f58SPeter Brune .  changed_y - indicates that the line search changed y
2358e557f58SPeter Brune .  changed_w - indicates that the line search changed w
23686d74e61SPeter Brune 
23786d74e61SPeter Brune    Level: intermediate
23886d74e61SPeter Brune 
23986d74e61SPeter Brune .keywords: set, linesearch, post-check
24086d74e61SPeter Brune 
241f1c6b773SPeter Brune .seealso: SNESLineSearchSetPreCheck()
24286d74e61SPeter Brune @*/
243f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetPostCheck(SNESLineSearch linesearch, SNESLineSearchPostCheckFunc func,void *ctx)
24486d74e61SPeter Brune {
24586d74e61SPeter Brune   PetscFunctionBegin;
246f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
24786d74e61SPeter Brune   if (func) linesearch->ops->postcheckstep = func;
24886d74e61SPeter Brune   if (ctx) linesearch->postcheckctx = ctx;
24986d74e61SPeter Brune   PetscFunctionReturn(0);
25086d74e61SPeter Brune }
25186d74e61SPeter Brune 
25286d74e61SPeter Brune 
25386d74e61SPeter Brune #undef __FUNCT__
254f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetPostCheck"
25586d74e61SPeter Brune /*@C
256f1c6b773SPeter Brune    SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine.
25786d74e61SPeter Brune 
25886d74e61SPeter Brune    Input Parameters:
259f1c6b773SPeter Brune .  linesearch - the SNESLineSearch context
26086d74e61SPeter Brune 
26186d74e61SPeter Brune    Output Parameters:
26286d74e61SPeter Brune +  func       - [optional] function evaluation routine
26386d74e61SPeter Brune -  ctx        - [optional] user-defined context for private data for the
26486d74e61SPeter Brune                 function evaluation routine (may be PETSC_NULL)
26586d74e61SPeter Brune 
26686d74e61SPeter Brune    Level: intermediate
26786d74e61SPeter Brune 
26886d74e61SPeter Brune .keywords: get, linesearch, post-check
26986d74e61SPeter Brune 
270f1c6b773SPeter Brune .seealso: SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck()
27186d74e61SPeter Brune @*/
272f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetPostCheck(SNESLineSearch linesearch, SNESLineSearchPostCheckFunc *func,void **ctx)
27386d74e61SPeter Brune {
2749bd66eb0SPeter Brune   PetscFunctionBegin;
275f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
27686d74e61SPeter Brune   if (func) *func = linesearch->ops->postcheckstep;
27786d74e61SPeter Brune   if (ctx) *ctx = linesearch->postcheckctx;
27886d74e61SPeter Brune   PetscFunctionReturn(0);
27986d74e61SPeter Brune }
28086d74e61SPeter Brune 
28186d74e61SPeter Brune 
282bf7f4e0aSPeter Brune #undef __FUNCT__
283f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPreCheck"
284f40b411bSPeter Brune /*@
285f1c6b773SPeter Brune    SNESLineSearchPreCheck - Prepares the line search for being applied.
286f40b411bSPeter Brune 
287cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
288f40b411bSPeter Brune 
289f40b411bSPeter Brune    Input Parameters:
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);
310*38bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed,4);
311bf7f4e0aSPeter Brune   }
312bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
313bf7f4e0aSPeter Brune }
314bf7f4e0aSPeter Brune 
315bf7f4e0aSPeter Brune #undef __FUNCT__
316f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPostCheck"
317f40b411bSPeter Brune /*@
318f1c6b773SPeter Brune    SNESLineSearchPostCheck - Prepares the line search for being applied.
319f40b411bSPeter Brune 
320cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
321f40b411bSPeter Brune 
322f40b411bSPeter Brune    Input Parameters:
3237b1df9c1SPeter Brune +  linesearch - The linesearch context
3247b1df9c1SPeter Brune .  X - The last solution
3257b1df9c1SPeter Brune .  Y - The step direction
3267b1df9c1SPeter Brune -  W - The updated solution, W = X + lambda*Y for some lambda
327f40b411bSPeter Brune 
328f40b411bSPeter Brune    Output Parameters:
32978bcb3b5SPeter Brune +  changed_Y - Indicator if the direction Y has been changed.
33078bcb3b5SPeter Brune -  changed_W - Indicator if the new candidate solution W has been changed.
331f40b411bSPeter Brune 
332f40b411bSPeter Brune    Level: Intermediate
333f40b411bSPeter Brune 
334f1c6b773SPeter Brune .keywords: SNESLineSearch, Create
335f40b411bSPeter Brune 
336f1c6b773SPeter Brune .seealso: SNESLineSearchPreCheck()
337f40b411bSPeter Brune @*/
3387b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch,Vec X,Vec Y,Vec W,PetscBool *changed_Y,PetscBool *changed_W)
339bf7f4e0aSPeter Brune {
340bf7f4e0aSPeter Brune   PetscErrorCode ierr;
341bf7f4e0aSPeter Brune   PetscFunctionBegin;
342bf7f4e0aSPeter Brune   *changed_Y = PETSC_FALSE;
343bf7f4e0aSPeter Brune   *changed_W = PETSC_FALSE;
344bf7f4e0aSPeter Brune   if (linesearch->ops->postcheckstep) {
3457b1df9c1SPeter Brune     ierr = (*linesearch->ops->postcheckstep)(linesearch,X,Y,W,changed_Y,changed_W,linesearch->postcheckctx);CHKERRQ(ierr);
346*38bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed_Y,5);
347*38bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed_W,6);
34886d74e61SPeter Brune   }
34986d74e61SPeter Brune   PetscFunctionReturn(0);
35086d74e61SPeter Brune }
35186d74e61SPeter Brune 
35286d74e61SPeter Brune 
35386d74e61SPeter Brune #undef __FUNCT__
354f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPreCheckPicard"
35586d74e61SPeter Brune /*@C
35686d74e61SPeter Brune    SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration
35786d74e61SPeter Brune 
358cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
35986d74e61SPeter Brune 
36086d74e61SPeter Brune    Input Arguments:
36186d74e61SPeter Brune +  linesearch - linesearch context
36286d74e61SPeter Brune .  X - base state for this step
36386d74e61SPeter Brune .  Y - initial correction
36486d74e61SPeter Brune 
36586d74e61SPeter Brune    Output Arguments:
36686d74e61SPeter Brune +  Y - correction, possibly modified
36786d74e61SPeter Brune -  changed - flag indicating that Y was modified
36886d74e61SPeter Brune 
36986d74e61SPeter Brune    Options Database Key:
370cd7522eaSPeter Brune +  -snes_linesearch_precheck_picard - activate this routine
371cd7522eaSPeter Brune -  -snes_linesearch_precheck_picard_angle - angle
37286d74e61SPeter Brune 
37386d74e61SPeter Brune    Level: advanced
37486d74e61SPeter Brune 
37586d74e61SPeter Brune    Notes:
37686d74e61SPeter Brune    This function should be passed to SNESLineSearchSetPreCheck()
37786d74e61SPeter Brune 
37886d74e61SPeter Brune    The justification for this method involves the linear convergence of a Picard iteration
37986d74e61SPeter Brune    so the Picard linearization should be provided in place of the "Jacobian". This correction
38086d74e61SPeter Brune    is generally not useful when using a Newton linearization.
38186d74e61SPeter Brune 
38286d74e61SPeter Brune    Reference:
38386d74e61SPeter Brune    Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology.
38486d74e61SPeter Brune 
38586d74e61SPeter Brune .seealso: SNESLineSearchSetPreCheck()
38686d74e61SPeter Brune @*/
387f1c6b773SPeter Brune PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed,void *ctx)
38886d74e61SPeter Brune {
38986d74e61SPeter Brune   PetscErrorCode ierr;
39086d74e61SPeter Brune   PetscReal      angle = *(PetscReal*)linesearch->precheckctx;
39186d74e61SPeter Brune   Vec            Ylast;
39286d74e61SPeter Brune   PetscScalar    dot;
393c87759e9SPeter Brune 
39486d74e61SPeter Brune   PetscInt       iter;
39586d74e61SPeter Brune   PetscReal      ynorm,ylastnorm,theta,angle_radians;
39686d74e61SPeter Brune   SNES           snes;
39786d74e61SPeter Brune 
39886d74e61SPeter Brune   PetscFunctionBegin;
399f1c6b773SPeter Brune   ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
40086d74e61SPeter Brune   ierr = PetscObjectQuery((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject*)&Ylast);CHKERRQ(ierr);
40186d74e61SPeter Brune   if (!Ylast) {
40286d74e61SPeter Brune     ierr = VecDuplicate(Y,&Ylast);CHKERRQ(ierr);
40386d74e61SPeter Brune     ierr = PetscObjectCompose((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject)Ylast);CHKERRQ(ierr);
40486d74e61SPeter Brune     ierr = PetscObjectDereference((PetscObject)Ylast);CHKERRQ(ierr);
40586d74e61SPeter Brune   }
40686d74e61SPeter Brune   ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr);
40786d74e61SPeter Brune   if (iter < 2) {
40886d74e61SPeter Brune     ierr = VecCopy(Y,Ylast);CHKERRQ(ierr);
40986d74e61SPeter Brune     *changed = PETSC_FALSE;
41086d74e61SPeter Brune     PetscFunctionReturn(0);
41186d74e61SPeter Brune   }
41286d74e61SPeter Brune 
41386d74e61SPeter Brune   ierr = VecDot(Y,Ylast,&dot);CHKERRQ(ierr);
41486d74e61SPeter Brune   ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr);
41586d74e61SPeter Brune   ierr = VecNorm(Ylast,NORM_2,&ylastnorm);CHKERRQ(ierr);
41686d74e61SPeter Brune   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
41786d74e61SPeter Brune   theta = acos((double)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0));
41886d74e61SPeter Brune   angle_radians = angle * PETSC_PI / 180.;
41986d74e61SPeter Brune   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
42086d74e61SPeter Brune     /* Modify the step Y */
42186d74e61SPeter Brune     PetscReal alpha,ydiffnorm;
42286d74e61SPeter Brune     ierr = VecAXPY(Ylast,-1.0,Y);CHKERRQ(ierr);
42386d74e61SPeter Brune     ierr = VecNorm(Ylast,NORM_2,&ydiffnorm);CHKERRQ(ierr);
42486d74e61SPeter Brune     alpha = ylastnorm / ydiffnorm;
42586d74e61SPeter Brune     ierr = VecCopy(Y,Ylast);CHKERRQ(ierr);
42686d74e61SPeter Brune     ierr = VecScale(Y,alpha);CHKERRQ(ierr);
42786d74e61SPeter 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);
42886d74e61SPeter Brune   } else {
42986d74e61SPeter Brune     ierr = PetscInfo2(snes,"Angle %G degrees exceeds threshold %G, no correction applied\n",theta*180./PETSC_PI,angle);CHKERRQ(ierr);
43086d74e61SPeter Brune     ierr = VecCopy(Y,Ylast);CHKERRQ(ierr);
43186d74e61SPeter Brune     *changed = PETSC_FALSE;
432bf7f4e0aSPeter Brune   }
433bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
434bf7f4e0aSPeter Brune }
435bf7f4e0aSPeter Brune 
436bf7f4e0aSPeter Brune #undef __FUNCT__
437f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchApply"
438f40b411bSPeter Brune /*@
439cd7522eaSPeter Brune    SNESLineSearchApply - Computes the line-search update.
440f40b411bSPeter Brune 
441f1c6b773SPeter Brune    Collective on SNESLineSearch
442f40b411bSPeter Brune 
443f40b411bSPeter Brune    Input Parameters:
4448e557f58SPeter Brune +  linesearch - The linesearch context
4458e557f58SPeter Brune .  X - The current solution
4468e557f58SPeter Brune .  F - The current function
4478e557f58SPeter Brune .  fnorm - The current norm
4488e557f58SPeter Brune -  Y - The search direction
449f40b411bSPeter Brune 
450f40b411bSPeter Brune    Output Parameters:
4518e557f58SPeter Brune +  X - The new solution
4528e557f58SPeter Brune .  F - The new function
4538e557f58SPeter Brune -  fnorm - The new function norm
454f40b411bSPeter Brune 
455cd7522eaSPeter Brune    Options Database Keys:
456cd7522eaSPeter Brune + -snes_linesearch_type - The Line search method
457cd7522eaSPeter Brune . -snes_linesearch_monitor - Print progress of line searches
458cd7522eaSPeter Brune . -snes_linesearch_damping - The linesearch damping parameter.
459cd7522eaSPeter Brune . -snes_linesearch_norms   - Turn on/off the linesearch norms
460cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess.
461cd7522eaSPeter Brune - -snes_linesearch_max_it - The number of iterations for iterative line searches.
462cd7522eaSPeter Brune 
463cd7522eaSPeter Brune    Notes:
464cd7522eaSPeter Brune    This is typically called from within a SNESSolve() implementation in order to
465cd7522eaSPeter Brune    help with convergence of the nonlinear method.  Various SNES types use line searches
466cd7522eaSPeter Brune    in different ways, but the overarching theme is that a line search is used to determine
467cd7522eaSPeter Brune    an optimal damping parameter of a step at each iteration of the method.  Each
468cd7522eaSPeter Brune    application of the line search may invoke SNESComputeFunction several times, and
469cd7522eaSPeter Brune    therefore may be fairly expensive.
470cd7522eaSPeter Brune 
471f40b411bSPeter Brune    Level: Intermediate
472f40b411bSPeter Brune 
473f1c6b773SPeter Brune .keywords: SNESLineSearch, Create
474f40b411bSPeter Brune 
475cd7522eaSPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchPreCheck(), SNESLineSearchPostCheck(), SNESSolve(), SNESComputeFunction()
476f40b411bSPeter Brune @*/
477f1c6b773SPeter Brune PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y) {
478bf7f4e0aSPeter Brune   PetscErrorCode ierr;
479bf7f4e0aSPeter Brune   PetscFunctionBegin;
480bf7f4e0aSPeter Brune 
481bf7f4e0aSPeter Brune   /* check the pointers */
482f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
483bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
484bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
485bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
486bf7f4e0aSPeter Brune 
487bf7f4e0aSPeter Brune   linesearch->success = PETSC_TRUE;
488bf7f4e0aSPeter Brune 
489bf7f4e0aSPeter Brune   linesearch->vec_sol = X;
490bf7f4e0aSPeter Brune   linesearch->vec_update = Y;
491bf7f4e0aSPeter Brune   linesearch->vec_func = F;
492bf7f4e0aSPeter Brune 
493f1c6b773SPeter Brune   ierr = SNESLineSearchSetUp(linesearch);CHKERRQ(ierr);
494bf7f4e0aSPeter Brune 
495bf7f4e0aSPeter Brune   if (!linesearch->keeplambda)
496bf7f4e0aSPeter Brune     linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */
497bf7f4e0aSPeter Brune 
498bf7f4e0aSPeter Brune   if (fnorm) {
499bf7f4e0aSPeter Brune     linesearch->fnorm = *fnorm;
500bf7f4e0aSPeter Brune   } else {
501bf7f4e0aSPeter Brune     ierr = VecNorm(F, NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
502bf7f4e0aSPeter Brune   }
503bf7f4e0aSPeter Brune 
504f1c6b773SPeter Brune   ierr = PetscLogEventBegin(SNESLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
505bf7f4e0aSPeter Brune 
506bf7f4e0aSPeter Brune   ierr = (*linesearch->ops->apply)(linesearch);CHKERRQ(ierr);
507bf7f4e0aSPeter Brune 
508f1c6b773SPeter Brune   ierr = PetscLogEventEnd(SNESLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
509bf7f4e0aSPeter Brune 
510bf7f4e0aSPeter Brune   if (fnorm)
511bf7f4e0aSPeter Brune     *fnorm = linesearch->fnorm;
512bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
513bf7f4e0aSPeter Brune }
514bf7f4e0aSPeter Brune 
515bf7f4e0aSPeter Brune #undef __FUNCT__
516f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchDestroy"
517f40b411bSPeter Brune /*@
518f1c6b773SPeter Brune    SNESLineSearchDestroy - Destroys the line search instance.
519f40b411bSPeter Brune 
520f1c6b773SPeter Brune    Collective on SNESLineSearch
521f40b411bSPeter Brune 
522f40b411bSPeter Brune    Input Parameters:
5238e557f58SPeter Brune .  linesearch - The linesearch context
524f40b411bSPeter Brune 
525f40b411bSPeter Brune    Level: Intermediate
526f40b411bSPeter Brune 
52778bcb3b5SPeter Brune .keywords: SNESLineSearch, Destroy
528f40b411bSPeter Brune 
529cd7522eaSPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchReset(), SNESDestroy()
530f40b411bSPeter Brune @*/
531f1c6b773SPeter Brune PetscErrorCode SNESLineSearchDestroy(SNESLineSearch * linesearch) {
532bf7f4e0aSPeter Brune   PetscErrorCode ierr;
533bf7f4e0aSPeter Brune   PetscFunctionBegin;
534bf7f4e0aSPeter Brune   if (!*linesearch) PetscFunctionReturn(0);
535f1c6b773SPeter Brune   PetscValidHeaderSpecific((*linesearch),SNESLINESEARCH_CLASSID,1);
536bf7f4e0aSPeter Brune   if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = 0; PetscFunctionReturn(0);}
537bf7f4e0aSPeter Brune   ierr = PetscObjectDepublish((*linesearch));CHKERRQ(ierr);
538f1c6b773SPeter Brune   ierr = SNESLineSearchReset(*linesearch);
539bf7f4e0aSPeter Brune   if ((*linesearch)->ops->destroy) {
540bf7f4e0aSPeter Brune     (*linesearch)->ops->destroy(*linesearch);
541bf7f4e0aSPeter Brune   }
542bf7f4e0aSPeter Brune   ierr = PetscViewerDestroy(&(*linesearch)->monitor);CHKERRQ(ierr);
543e7058c64SPeter Brune   ierr = PetscHeaderDestroy(linesearch);CHKERRQ(ierr);
544bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
545bf7f4e0aSPeter Brune }
546bf7f4e0aSPeter Brune 
547bf7f4e0aSPeter Brune #undef __FUNCT__
548f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetMonitor"
549f40b411bSPeter Brune /*@
550cd7522eaSPeter Brune    SNESLineSearchSetMonitor - Turns on/off printing useful information and debugging output about the line search.
551bf7f4e0aSPeter Brune 
552bf7f4e0aSPeter Brune    Input Parameters:
553bf7f4e0aSPeter Brune +  snes - nonlinear context obtained from SNESCreate()
554bf7f4e0aSPeter Brune -  flg - PETSC_TRUE to monitor the line search
555bf7f4e0aSPeter Brune 
556bf7f4e0aSPeter Brune    Logically Collective on SNES
557bf7f4e0aSPeter Brune 
558bf7f4e0aSPeter Brune    Options Database:
5598e557f58SPeter Brune .   -snes_linesearch_monitor - enables the monitor
560bf7f4e0aSPeter Brune 
561bf7f4e0aSPeter Brune    Level: intermediate
562bf7f4e0aSPeter Brune 
563bf7f4e0aSPeter Brune 
564cd7522eaSPeter Brune .seealso: SNESLineSearchGetMonitor(), PetscViewer
565bf7f4e0aSPeter Brune @*/
566f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetMonitor(SNESLineSearch linesearch, PetscBool flg)
567bf7f4e0aSPeter Brune {
568bf7f4e0aSPeter Brune 
569bf7f4e0aSPeter Brune   PetscErrorCode ierr;
570bf7f4e0aSPeter Brune   PetscFunctionBegin;
571bf7f4e0aSPeter Brune   if (flg && !linesearch->monitor) {
572bf7f4e0aSPeter Brune     ierr = PetscViewerASCIIOpen(((PetscObject)linesearch)->comm,"stdout",&linesearch->monitor);CHKERRQ(ierr);
573bf7f4e0aSPeter Brune   } else if (!flg && linesearch->monitor) {
574bf7f4e0aSPeter Brune     ierr = PetscViewerDestroy(&linesearch->monitor);CHKERRQ(ierr);
575bf7f4e0aSPeter Brune   }
576bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
577bf7f4e0aSPeter Brune }
578bf7f4e0aSPeter Brune 
579bf7f4e0aSPeter Brune #undef __FUNCT__
580f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetMonitor"
581f40b411bSPeter Brune /*@
582cd7522eaSPeter Brune    SNESLineSearchGetMonitor - Gets the PetscViewer instance for the line search monitor.
5836a388c36SPeter Brune 
584f40b411bSPeter Brune    Input Parameters:
5858e557f58SPeter Brune .  linesearch - linesearch context
586f40b411bSPeter Brune 
587f40b411bSPeter Brune    Input Parameters:
5888e557f58SPeter Brune .  monitor - monitor context
589f40b411bSPeter Brune 
590f40b411bSPeter Brune    Logically Collective on SNES
591f40b411bSPeter Brune 
592f40b411bSPeter Brune 
593f40b411bSPeter Brune    Options Database Keys:
5948e557f58SPeter Brune .   -snes_linesearch_monitor - enables the monitor
595f40b411bSPeter Brune 
596f40b411bSPeter Brune    Level: intermediate
597f40b411bSPeter Brune 
598f40b411bSPeter Brune 
599cd7522eaSPeter Brune .seealso: SNESLineSearchSetMonitor(), PetscViewer
600f40b411bSPeter Brune @*/
601f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetMonitor(SNESLineSearch linesearch, PetscViewer *monitor)
6026a388c36SPeter Brune {
6036a388c36SPeter Brune 
6046a388c36SPeter Brune   PetscFunctionBegin;
605f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
6066a388c36SPeter Brune   if (monitor) {
6076a388c36SPeter Brune     PetscValidPointer(monitor, 2);
6086a388c36SPeter Brune     *monitor = linesearch->monitor;
6096a388c36SPeter Brune   }
6106a388c36SPeter Brune   PetscFunctionReturn(0);
6116a388c36SPeter Brune }
6126a388c36SPeter Brune 
6136a388c36SPeter Brune #undef __FUNCT__
614f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetFromOptions"
615f40b411bSPeter Brune /*@
616f1c6b773SPeter Brune    SNESLineSearchSetFromOptions - Sets options for the line search
617f40b411bSPeter Brune 
618f40b411bSPeter Brune    Input Parameters:
6198e557f58SPeter Brune .  linesearch - linesearch context
620f40b411bSPeter Brune 
621f40b411bSPeter Brune    Options Database Keys:
622cd7522eaSPeter Brune + -snes_linesearch_type - The Line search method
62371eef1aeSPeter Brune . -snes_linesearch_minlambda - The minimum step length
6241a9b3a06SPeter Brune . -snes_linesearch_maxstep - The maximum step size
6251a9b3a06SPeter Brune . -snes_linesearch_rtol - Relative tolerance for iterative line searches
6261a9b3a06SPeter Brune . -snes_linesearch_atol - Absolute tolerance for iterative line searches
6271a9b3a06SPeter Brune . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches
6281a9b3a06SPeter Brune . -snes_linesearch_max_it - The number of iterations for iterative line searches
629cd7522eaSPeter Brune . -snes_linesearch_monitor - Print progress of line searches
6308e557f58SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter
631cd7522eaSPeter Brune . -snes_linesearch_norms   - Turn on/off the linesearch norms
632cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess.
63359405d5eSPeter Brune . -snes_linesearch_order - The polynomial order of approximation used in the linesearch
6341a9b3a06SPeter Brune . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method
6351a9b3a06SPeter Brune - -snes_linesearch_precheck_picard_angle - Angle used in picard precheck method
636f40b411bSPeter Brune 
637f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
638f40b411bSPeter Brune 
639f40b411bSPeter Brune    Level: intermediate
640f40b411bSPeter Brune 
641f40b411bSPeter Brune 
642f1c6b773SPeter Brune .seealso: SNESLineSearchCreate()
643f40b411bSPeter Brune @*/
644f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch) {
645bf7f4e0aSPeter Brune   PetscErrorCode ierr;
6461a4f838cSPeter Brune   const char     *deft = SNESLINESEARCHBASIC;
647bf7f4e0aSPeter Brune   char           type[256];
648bf7f4e0aSPeter Brune   PetscBool      flg, set;
649bf7f4e0aSPeter Brune   PetscFunctionBegin;
650f1c6b773SPeter Brune   if (!SNESLineSearchRegisterAllCalled) {ierr = SNESLineSearchRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
651bf7f4e0aSPeter Brune 
652bf7f4e0aSPeter Brune   ierr = PetscObjectOptionsBegin((PetscObject)linesearch);CHKERRQ(ierr);
653bf7f4e0aSPeter Brune   if (((PetscObject)linesearch)->type_name) {
654bf7f4e0aSPeter Brune     deft = ((PetscObject)linesearch)->type_name;
655bf7f4e0aSPeter Brune   }
6567a35526eSPeter Brune   ierr = PetscOptionsList("-snes_linesearch_type","Line-search method","SNESLineSearchSetType",SNESLineSearchList,deft,type,256,&flg);CHKERRQ(ierr);
657bf7f4e0aSPeter Brune   if (flg) {
658f1c6b773SPeter Brune     ierr = SNESLineSearchSetType(linesearch,type);CHKERRQ(ierr);
659bf7f4e0aSPeter Brune   } else if (!((PetscObject)linesearch)->type_name) {
660f1c6b773SPeter Brune     ierr = SNESLineSearchSetType(linesearch,deft);CHKERRQ(ierr);
661bf7f4e0aSPeter Brune   }
662bf7f4e0aSPeter Brune   if (linesearch->ops->setfromoptions) {
663bf7f4e0aSPeter Brune     (*linesearch->ops->setfromoptions)(linesearch);CHKERRQ(ierr);
664bf7f4e0aSPeter Brune   }
665bf7f4e0aSPeter Brune 
6667a35526eSPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_monitor","Print progress of line searches","SNESSNESLineSearchSetMonitor",
667bf7f4e0aSPeter Brune                           linesearch->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
668f1c6b773SPeter Brune   if (set) {ierr = SNESLineSearchSetMonitor(linesearch,flg);CHKERRQ(ierr);}
669bf7f4e0aSPeter Brune 
6701a9b3a06SPeter Brune   /* tolerances */
67171eef1aeSPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_minlambda","Minimum step length","SNESLineSearchSetTolerances",linesearch->steptol,&linesearch->steptol,0);CHKERRQ(ierr);
6721a9b3a06SPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_maxstep","Maximum step size","SNESLineSearchSetTolerances",linesearch->maxstep,&linesearch->maxstep,0);CHKERRQ(ierr);
6731a9b3a06SPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_rtol","Relative tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->rtol,&linesearch->rtol,0);CHKERRQ(ierr);
6741a9b3a06SPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_atol","Absolute tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->atol,&linesearch->atol,0);CHKERRQ(ierr);
6751a9b3a06SPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_ltol","Change in lambda tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->ltol,&linesearch->ltol,0);CHKERRQ(ierr);
6768e557f58SPeter Brune   ierr = PetscOptionsInt("-snes_linesearch_max_it","Maximum iterations for iterative line searches","SNESLineSearchSetTolerances",linesearch->max_its,&linesearch->max_its,0);CHKERRQ(ierr);
6777a35526eSPeter Brune 
6781a9b3a06SPeter Brune   /* damping parameters */
6791a9b3a06SPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_damping","Line search damping and initial step guess","SNESLineSearchSetDamping",linesearch->damping,&linesearch->damping,0);CHKERRQ(ierr);
6801a9b3a06SPeter Brune 
6811a9b3a06SPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_keeplambda","Use previous lambda as damping","SNESLineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,0);CHKERRQ(ierr);
6821a9b3a06SPeter Brune 
6831a9b3a06SPeter Brune   /* precheck */
6847a35526eSPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_precheck_picard","Use a correction that sometimes improves convergence of Picard iteration","SNESLineSearchPreCheckPicard",flg,&flg,&set);CHKERRQ(ierr);
6857a35526eSPeter Brune   if (set) {
6867a35526eSPeter Brune     if (flg) {
6877a35526eSPeter Brune       linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */
6887a35526eSPeter Brune       ierr = PetscOptionsReal("-snes_linesearch_precheck_picard_angle","Maximum angle at which to activate the correction",
6897a35526eSPeter Brune                               "none",linesearch->precheck_picard_angle,&linesearch->precheck_picard_angle,PETSC_NULL);CHKERRQ(ierr);
6907a35526eSPeter Brune       ierr = SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&linesearch->precheck_picard_angle);CHKERRQ(ierr);
6917a35526eSPeter Brune     } else {
6927a35526eSPeter Brune       ierr = SNESLineSearchSetPreCheck(linesearch,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
6937a35526eSPeter Brune     }
6947a35526eSPeter Brune   }
695b000cd8dSPeter Brune   ierr = PetscOptionsInt("-snes_linesearch_order","Order of approximation used in the line search","SNESLineSearchSetOrder",linesearch->order,&linesearch->order,0);CHKERRQ(ierr);
6961a9b3a06SPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_norms","Compute final norms in line search","SNESLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,0);CHKERRQ(ierr);
6977a35526eSPeter Brune 
698bf7f4e0aSPeter Brune   ierr = PetscObjectProcessOptionsHandlers((PetscObject)linesearch);CHKERRQ(ierr);
699bf7f4e0aSPeter Brune   ierr = PetscOptionsEnd();CHKERRQ(ierr);
700bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
701bf7f4e0aSPeter Brune }
702bf7f4e0aSPeter Brune 
703bf7f4e0aSPeter Brune #undef __FUNCT__
704f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchView"
705f40b411bSPeter Brune /*@
706cd7522eaSPeter Brune    SNESLineSearchView - Prints useful information about the line search not
707cd7522eaSPeter Brune    related to an individual call.
708f40b411bSPeter Brune 
709f40b411bSPeter Brune    Input Parameters:
7108e557f58SPeter Brune .  linesearch - linesearch context
711f40b411bSPeter Brune 
712f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
713f40b411bSPeter Brune 
714f40b411bSPeter Brune    Level: intermediate
715f40b411bSPeter Brune 
716f1c6b773SPeter Brune .seealso: SNESLineSearchCreate()
717f40b411bSPeter Brune @*/
7187f1410a3SPeter Brune PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer) {
7197f1410a3SPeter Brune   PetscErrorCode ierr;
7207f1410a3SPeter Brune   PetscBool      iascii;
721bf7f4e0aSPeter Brune   PetscFunctionBegin;
7227f1410a3SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
7237f1410a3SPeter Brune   if (!viewer) {
7247f1410a3SPeter Brune     ierr = PetscViewerASCIIGetStdout(((PetscObject)linesearch)->comm,&viewer);CHKERRQ(ierr);
7257f1410a3SPeter Brune   }
7267f1410a3SPeter Brune   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
7277f1410a3SPeter Brune   PetscCheckSameComm(linesearch,1,viewer,2);
728f40b411bSPeter Brune 
729251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
7307f1410a3SPeter Brune   if (iascii) {
7317f1410a3SPeter Brune     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)linesearch,viewer,"SNESLineSearch Object");CHKERRQ(ierr);
7327f1410a3SPeter Brune     if (linesearch->ops->view) {
7337f1410a3SPeter Brune       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
7347f1410a3SPeter Brune       ierr = (*linesearch->ops->view)(linesearch,viewer);CHKERRQ(ierr);
7357f1410a3SPeter Brune       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
7367f1410a3SPeter Brune     }
737007b6e36SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  maxstep=%e, minlambda=%e\n", linesearch->maxstep,linesearch->steptol);CHKERRQ(ierr);
738007b6e36SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%e, absolute=%e, lambda=%e\n", linesearch->rtol,linesearch->atol,linesearch->ltol);CHKERRQ(ierr);
7397f1410a3SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D\n", linesearch->max_its);CHKERRQ(ierr);
7407f1410a3SPeter Brune     if (linesearch->ops->precheckstep) {
7417f1410a3SPeter Brune       if (linesearch->ops->precheckstep == SNESLineSearchPreCheckPicard) {
7427f1410a3SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  using precheck step to speed up Picard convergence\n", linesearch->max_its);CHKERRQ(ierr);
7437f1410a3SPeter Brune       } else {
7447f1410a3SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  using user-defined precheck step\n", linesearch->max_its);CHKERRQ(ierr);
7457f1410a3SPeter Brune       }
7467f1410a3SPeter Brune     }
7477f1410a3SPeter Brune     if (linesearch->ops->postcheckstep) {
7487f1410a3SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer,"  using user-defined postcheck step\n", linesearch->max_its);CHKERRQ(ierr);
7497f1410a3SPeter Brune     }
7507f1410a3SPeter Brune   }
751bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
752bf7f4e0aSPeter Brune }
753bf7f4e0aSPeter Brune 
754bf7f4e0aSPeter Brune #undef __FUNCT__
755f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetType"
756ea5d4fccSPeter Brune /*@C
757f1c6b773SPeter Brune    SNESLineSearchSetType - Sets the linesearch type
758f40b411bSPeter Brune 
759f40b411bSPeter Brune    Input Parameters:
7608e557f58SPeter Brune +  linesearch - linesearch context
761f40b411bSPeter Brune -  type - The type of line search to be used
762f40b411bSPeter Brune 
763cd7522eaSPeter Brune    Available Types:
764cd7522eaSPeter Brune +  basic - Simple damping line search.
7658e557f58SPeter Brune .  bt - Backtracking line search over the L2 norm of the function
7668e557f58SPeter Brune .  l2 - Secant line search over the L2 norm of the function
7678e557f58SPeter Brune .  cp - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x)
7688e557f58SPeter Brune -  shell - User provided SNESLineSearch implementation
769cd7522eaSPeter Brune 
770f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
771f40b411bSPeter Brune 
772f40b411bSPeter Brune    Level: intermediate
773f40b411bSPeter Brune 
774f40b411bSPeter Brune 
775f1c6b773SPeter Brune .seealso: SNESLineSearchCreate()
776f40b411bSPeter Brune @*/
777f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, const SNESLineSearchType type)
778bf7f4e0aSPeter Brune {
779bf7f4e0aSPeter Brune 
780f1c6b773SPeter Brune   PetscErrorCode ierr,(*r)(SNESLineSearch);
781bf7f4e0aSPeter Brune   PetscBool      match;
782bf7f4e0aSPeter Brune 
783bf7f4e0aSPeter Brune   PetscFunctionBegin;
784f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
785bf7f4e0aSPeter Brune   PetscValidCharPointer(type,2);
786bf7f4e0aSPeter Brune 
787251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)linesearch,type,&match);CHKERRQ(ierr);
788bf7f4e0aSPeter Brune   if (match) PetscFunctionReturn(0);
789bf7f4e0aSPeter Brune 
790f1c6b773SPeter Brune   ierr =  PetscFListFind(SNESLineSearchList,((PetscObject)linesearch)->comm,type,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr);
791bf7f4e0aSPeter Brune   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type);
792bf7f4e0aSPeter Brune   /* Destroy the previous private linesearch context */
793bf7f4e0aSPeter Brune   if (linesearch->ops->destroy) {
794bf7f4e0aSPeter Brune     ierr = (*(linesearch)->ops->destroy)(linesearch);CHKERRQ(ierr);
795bf7f4e0aSPeter Brune     linesearch->ops->destroy = PETSC_NULL;
796bf7f4e0aSPeter Brune   }
797f1c6b773SPeter Brune   /* Reinitialize function pointers in SNESLineSearchOps structure */
798bf7f4e0aSPeter Brune   linesearch->ops->apply          = 0;
799bf7f4e0aSPeter Brune   linesearch->ops->view           = 0;
800bf7f4e0aSPeter Brune   linesearch->ops->setfromoptions = 0;
801bf7f4e0aSPeter Brune   linesearch->ops->destroy        = 0;
802bf7f4e0aSPeter Brune 
803bf7f4e0aSPeter Brune   ierr = PetscObjectChangeTypeName((PetscObject)linesearch,type);CHKERRQ(ierr);
804bf7f4e0aSPeter Brune   ierr = (*r)(linesearch);CHKERRQ(ierr);
805bf7f4e0aSPeter Brune #if defined(PETSC_HAVE_AMS)
806bf7f4e0aSPeter Brune   if (PetscAMSPublishAll) {
807bf7f4e0aSPeter Brune     ierr = PetscObjectAMSPublish((PetscObject)linesearch);CHKERRQ(ierr);
808bf7f4e0aSPeter Brune   }
809bf7f4e0aSPeter Brune #endif
810bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
811bf7f4e0aSPeter Brune }
812bf7f4e0aSPeter Brune 
813bf7f4e0aSPeter Brune #undef __FUNCT__
814f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetSNES"
815f40b411bSPeter Brune /*@
81678bcb3b5SPeter Brune    SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation.
817f40b411bSPeter Brune 
818f40b411bSPeter Brune    Input Parameters:
8198e557f58SPeter Brune +  linesearch - linesearch context
820f40b411bSPeter Brune -  snes - The snes instance
821f40b411bSPeter Brune 
82278bcb3b5SPeter Brune    Level: developer
82378bcb3b5SPeter Brune 
82478bcb3b5SPeter Brune    Notes:
82578bcb3b5SPeter Brune    This happens automatically when the line search is gotten/created with
82678bcb3b5SPeter Brune    SNESGetSNESLineSearch().  This routine is therefore mainly called within SNES
82778bcb3b5SPeter Brune    implementations.
828f40b411bSPeter Brune 
8298141a3b9SPeter Brune    Level: developer
830f40b411bSPeter Brune 
831cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
832f40b411bSPeter Brune @*/
833f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes){
834bf7f4e0aSPeter Brune   PetscFunctionBegin;
835f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
836bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,2);
837bf7f4e0aSPeter Brune   linesearch->snes = snes;
838bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
839bf7f4e0aSPeter Brune }
840bf7f4e0aSPeter Brune 
841bf7f4e0aSPeter Brune #undef __FUNCT__
842f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetSNES"
843f40b411bSPeter Brune /*@
8448141a3b9SPeter Brune    SNESLineSearchGetSNES - Gets the SNES instance associated with the line search.
8458141a3b9SPeter Brune    Having an associated SNES is necessary because most line search implementations must be able to
8468141a3b9SPeter Brune    evaluate the function using SNESComputeFunction() for the associated SNES.  This routine
8478141a3b9SPeter Brune    is used in the line search implementations when one must get this associated SNES instance.
848f40b411bSPeter Brune 
849f40b411bSPeter Brune    Input Parameters:
8508e557f58SPeter Brune .  linesearch - linesearch context
851f40b411bSPeter Brune 
852f40b411bSPeter Brune    Output Parameters:
853f40b411bSPeter Brune .  snes - The snes instance
854f40b411bSPeter Brune 
8558141a3b9SPeter Brune    Level: developer
856f40b411bSPeter Brune 
857cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
858f40b411bSPeter Brune @*/
859f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes){
860bf7f4e0aSPeter Brune   PetscFunctionBegin;
861f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
8626a388c36SPeter Brune   PetscValidPointer(snes, 2);
863bf7f4e0aSPeter Brune   *snes = linesearch->snes;
864bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
865bf7f4e0aSPeter Brune }
866bf7f4e0aSPeter Brune 
8676a388c36SPeter Brune #undef __FUNCT__
868f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetLambda"
869f40b411bSPeter Brune /*@
870f1c6b773SPeter Brune    SNESLineSearchGetLambda - Gets the last linesearch steplength discovered.
871f40b411bSPeter Brune 
872f40b411bSPeter Brune    Input Parameters:
8738e557f58SPeter Brune .  linesearch - linesearch context
874f40b411bSPeter Brune 
875f40b411bSPeter Brune    Output Parameters:
876cd7522eaSPeter Brune .  lambda - The last steplength computed during SNESLineSearchApply()
877f40b411bSPeter Brune 
87878bcb3b5SPeter Brune    Level: advanced
87978bcb3b5SPeter Brune 
8808e557f58SPeter Brune    Notes:
8818e557f58SPeter Brune    This is useful in methods where the solver is ill-scaled and
88278bcb3b5SPeter Brune    requires some adaptive notion of the difference in scale between the
88378bcb3b5SPeter Brune    solution and the function.  For instance, SNESQN may be scaled by the
88478bcb3b5SPeter Brune    line search lambda using the argument -snes_qn_scaling ls.
88578bcb3b5SPeter Brune 
886f40b411bSPeter Brune 
887cd7522eaSPeter Brune .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetDamping(), SNESLineSearchApply()
888f40b411bSPeter Brune @*/
889f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda)
8906a388c36SPeter Brune {
8916a388c36SPeter Brune   PetscFunctionBegin;
892f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
8936a388c36SPeter Brune   PetscValidPointer(lambda, 2);
8946a388c36SPeter Brune   *lambda = linesearch->lambda;
8956a388c36SPeter Brune   PetscFunctionReturn(0);
8966a388c36SPeter Brune }
8976a388c36SPeter Brune 
8986a388c36SPeter Brune #undef __FUNCT__
899f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetLambda"
900f40b411bSPeter Brune /*@
901f1c6b773SPeter Brune    SNESLineSearchSetLambda - Sets the linesearch steplength.
902f40b411bSPeter Brune 
903f40b411bSPeter Brune    Input Parameters:
9048e557f58SPeter Brune +  linesearch - linesearch context
905f40b411bSPeter Brune -  lambda - The last steplength.
906f40b411bSPeter Brune 
907cd7522eaSPeter Brune    Notes:
908cd7522eaSPeter Brune    This routine is typically used within implementations of SNESLineSearchApply
909cd7522eaSPeter Brune    to set the final steplength.  This routine (and SNESLineSearchGetLambda()) were
910cd7522eaSPeter Brune    added in order to facilitate Quasi-Newton methods that use the previous steplength
911cd7522eaSPeter Brune    as an inner scaling parameter.
912cd7522eaSPeter Brune 
91378bcb3b5SPeter Brune    Level: advanced
914f40b411bSPeter Brune 
915f1c6b773SPeter Brune .seealso: SNESLineSearchGetLambda()
916f40b411bSPeter Brune @*/
917f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
9186a388c36SPeter Brune {
9196a388c36SPeter Brune   PetscFunctionBegin;
920f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
9216a388c36SPeter Brune   linesearch->lambda = lambda;
9226a388c36SPeter Brune   PetscFunctionReturn(0);
9236a388c36SPeter Brune }
9246a388c36SPeter Brune 
9256a388c36SPeter Brune #undef  __FUNCT__
926f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetTolerances"
927f40b411bSPeter Brune /*@
92878bcb3b5SPeter Brune    SNESLineSearchGetTolerances - Gets the tolerances for the method.  These include
92978bcb3b5SPeter Brune    tolerances for the relative and absolute change in the function norm, the change
93078bcb3b5SPeter Brune    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
93178bcb3b5SPeter Brune    and the maximum number of iterations the line search procedure may take.
932f40b411bSPeter Brune 
933f40b411bSPeter Brune    Input Parameters:
9348e557f58SPeter Brune .  linesearch - linesearch context
935f40b411bSPeter Brune 
936f40b411bSPeter Brune    Output Parameters:
937516fe3c3SPeter Brune +  steptol - The minimum steplength
9386cc8e53bSPeter Brune .  maxstep - The maximum steplength
939516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
940516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
941516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
942516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
943f40b411bSPeter Brune 
94478bcb3b5SPeter Brune    Level: intermediate
94578bcb3b5SPeter Brune 
94678bcb3b5SPeter Brune    Notes:
94778bcb3b5SPeter Brune    Different line searches may implement these parameters slightly differently as
94878bcb3b5SPeter Brune    the method requires.
949516fe3c3SPeter Brune 
950f1c6b773SPeter Brune .seealso: SNESLineSearchSetTolerances()
951f40b411bSPeter Brune @*/
952f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
9536a388c36SPeter Brune {
9546a388c36SPeter Brune   PetscFunctionBegin;
955f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
956516fe3c3SPeter Brune   if (steptol) {
9576a388c36SPeter Brune     PetscValidPointer(steptol, 2);
9586a388c36SPeter Brune     *steptol = linesearch->steptol;
959516fe3c3SPeter Brune   }
960516fe3c3SPeter Brune   if (maxstep) {
961516fe3c3SPeter Brune     PetscValidPointer(maxstep, 3);
962516fe3c3SPeter Brune     *maxstep = linesearch->maxstep;
963516fe3c3SPeter Brune   }
964516fe3c3SPeter Brune   if (rtol) {
965516fe3c3SPeter Brune     PetscValidPointer(rtol, 4);
966516fe3c3SPeter Brune     *rtol = linesearch->rtol;
967516fe3c3SPeter Brune   }
968516fe3c3SPeter Brune   if (atol) {
969516fe3c3SPeter Brune     PetscValidPointer(atol, 5);
970516fe3c3SPeter Brune     *atol = linesearch->atol;
971516fe3c3SPeter Brune   }
972516fe3c3SPeter Brune   if (ltol) {
973516fe3c3SPeter Brune     PetscValidPointer(ltol, 6);
974516fe3c3SPeter Brune     *ltol = linesearch->ltol;
975516fe3c3SPeter Brune   }
976516fe3c3SPeter Brune   if (max_its) {
977516fe3c3SPeter Brune     PetscValidPointer(max_its, 7);
978516fe3c3SPeter Brune     *max_its = linesearch->max_its;
979516fe3c3SPeter Brune   }
9806a388c36SPeter Brune   PetscFunctionReturn(0);
9816a388c36SPeter Brune }
9826a388c36SPeter Brune 
9836a388c36SPeter Brune #undef  __FUNCT__
984f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetTolerances"
985f40b411bSPeter Brune /*@
98678bcb3b5SPeter Brune    SNESLineSearchSetTolerances -  Gets the tolerances for the method.  These include
98778bcb3b5SPeter Brune    tolerances for the relative and absolute change in the function norm, the change
98878bcb3b5SPeter Brune    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
98978bcb3b5SPeter Brune    and the maximum number of iterations the line search procedure may take.
990f40b411bSPeter Brune 
991f40b411bSPeter Brune    Input Parameters:
9928e557f58SPeter Brune +  linesearch - linesearch context
993516fe3c3SPeter Brune .  steptol - The minimum steplength
9946cc8e53bSPeter Brune .  maxstep - The maximum steplength
995516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
996516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
997516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
998516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
999f40b411bSPeter Brune 
100078bcb3b5SPeter Brune    Notes:
100178bcb3b5SPeter Brune    The user may choose to not set any of the tolerances in the method using PETSC_DEFAULT in
100278bcb3b5SPeter Brune    place of an argument.
1003f40b411bSPeter Brune 
100478bcb3b5SPeter Brune    Level: intermediate
1005516fe3c3SPeter Brune 
1006f1c6b773SPeter Brune .seealso: SNESLineSearchGetTolerances()
1007f40b411bSPeter Brune @*/
1008f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its)
10096a388c36SPeter Brune {
10106a388c36SPeter Brune   PetscFunctionBegin;
1011f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1012d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,steptol,2);
1013d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,maxstep,3);
1014d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,rtol,4);
1015d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,atol,5);
1016d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,ltol,6);
1017d3952184SSatish Balay   PetscValidLogicalCollectiveInt(linesearch,max_its,7);
1018d3952184SSatish Balay 
1019d3952184SSatish Balay   if ( steptol!= PETSC_DEFAULT) {
1020d3952184SSatish Balay     if ( steptol < 0.0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Minimum step length %G must be non-negative",steptol);
10216a388c36SPeter Brune     linesearch->steptol = steptol;
1022d3952184SSatish Balay   }
1023d3952184SSatish Balay 
1024d3952184SSatish Balay   if ( maxstep!= PETSC_DEFAULT) {
1025d3952184SSatish Balay     if ( maxstep < 0.0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Maximum step length %G must be non-negative",maxstep);
1026516fe3c3SPeter Brune     linesearch->maxstep = maxstep;
1027d3952184SSatish Balay   }
1028d3952184SSatish Balay 
1029d3952184SSatish Balay   if (rtol != PETSC_DEFAULT) {
1030d3952184SSatish 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);
1031516fe3c3SPeter Brune     linesearch->rtol = rtol;
1032d3952184SSatish Balay   }
1033d3952184SSatish Balay 
1034d3952184SSatish Balay   if (atol != PETSC_DEFAULT) {
1035d3952184SSatish Balay     if (atol < 0.0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %G must be non-negative",atol);
1036516fe3c3SPeter Brune     linesearch->atol = atol;
1037d3952184SSatish Balay   }
1038d3952184SSatish Balay 
1039d3952184SSatish Balay   if (ltol != PETSC_DEFAULT) {
1040d3952184SSatish Balay     if (ltol < 0.0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Labmda tolerance %G must be non-negative",ltol);
1041516fe3c3SPeter Brune   linesearch->ltol = ltol;
1042d3952184SSatish Balay   }
1043d3952184SSatish Balay 
1044d3952184SSatish Balay   if (max_its != PETSC_DEFAULT) {
1045d3952184SSatish Balay     if (max_its < 0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",max_its);
1046516fe3c3SPeter Brune     linesearch->max_its = max_its;
1047d3952184SSatish Balay   }
1048d3952184SSatish Balay 
10496a388c36SPeter Brune   PetscFunctionReturn(0);
10506a388c36SPeter Brune }
10516a388c36SPeter Brune 
1052516fe3c3SPeter Brune 
10536a388c36SPeter Brune #undef __FUNCT__
1054f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetDamping"
1055f40b411bSPeter Brune /*@
1056f1c6b773SPeter Brune    SNESLineSearchGetDamping - Gets the line search damping parameter.
1057f40b411bSPeter Brune 
1058f40b411bSPeter Brune    Input Parameters:
10598e557f58SPeter Brune .  linesearch - linesearch context
1060f40b411bSPeter Brune 
1061f40b411bSPeter Brune    Output Parameters:
10628e557f58SPeter Brune .  damping - The damping parameter
1063f40b411bSPeter Brune 
106478bcb3b5SPeter Brune    Level: advanced
1065f40b411bSPeter Brune 
106678bcb3b5SPeter Brune .seealso: SNESLineSearchGetStepTolerance(), SNESQN
1067f40b411bSPeter Brune @*/
1068f40b411bSPeter Brune 
1069f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping)
10706a388c36SPeter Brune {
10716a388c36SPeter Brune   PetscFunctionBegin;
1072f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
10736a388c36SPeter Brune   PetscValidPointer(damping, 2);
10746a388c36SPeter Brune   *damping = linesearch->damping;
10756a388c36SPeter Brune   PetscFunctionReturn(0);
10766a388c36SPeter Brune }
10776a388c36SPeter Brune 
10786a388c36SPeter Brune #undef __FUNCT__
1079f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetDamping"
1080f40b411bSPeter Brune /*@
1081f1c6b773SPeter Brune    SNESLineSearchSetDamping - Sets the line search damping paramter.
1082f40b411bSPeter Brune 
1083f40b411bSPeter Brune    Input Parameters:
108478bcb3b5SPeter Brune .  linesearch - linesearch context
108578bcb3b5SPeter Brune .  damping - The damping parameter
1086f40b411bSPeter Brune 
1087f40b411bSPeter Brune    Level: intermediate
1088f40b411bSPeter Brune 
1089cd7522eaSPeter Brune    Notes:
1090cd7522eaSPeter Brune    The basic line search merely takes the update step scaled by the damping parameter.
1091cd7522eaSPeter Brune    The use of the damping parameter in the l2 and cp line searches is much more subtle;
109278bcb3b5SPeter Brune    it is used as a starting point in calculating the secant step. However, the eventual
1093cd7522eaSPeter Brune    step may be of greater length than the damping parameter.  In the bt line search it is
1094cd7522eaSPeter Brune    used as the maximum possible step length, as the bt line search only backtracks.
1095cd7522eaSPeter Brune 
1096f1c6b773SPeter Brune .seealso: SNESLineSearchGetDamping()
1097f40b411bSPeter Brune @*/
1098f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping)
10996a388c36SPeter Brune {
11006a388c36SPeter Brune   PetscFunctionBegin;
1101f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
11026a388c36SPeter Brune   linesearch->damping = damping;
11036a388c36SPeter Brune   PetscFunctionReturn(0);
11046a388c36SPeter Brune }
11056a388c36SPeter Brune 
11066a388c36SPeter Brune #undef __FUNCT__
110759405d5eSPeter Brune #define __FUNCT__ "SNESLineSearchGetOrder"
110859405d5eSPeter Brune /*@
110959405d5eSPeter Brune    SNESLineSearchGetOrder - Gets the line search approximation order.
111059405d5eSPeter Brune 
111159405d5eSPeter Brune    Input Parameters:
111278bcb3b5SPeter Brune .  linesearch - linesearch context
111359405d5eSPeter Brune 
111459405d5eSPeter Brune    Output Parameters:
11158e557f58SPeter Brune .  order - The order
111659405d5eSPeter Brune 
111778bcb3b5SPeter Brune    Possible Values for order:
11188e557f58SPeter Brune +  SNES_LINESEARCH_LINEAR - linear method
11198e557f58SPeter Brune .  SNES_LINESEARCH_QUADRATIC - quadratic method
11208e557f58SPeter Brune -  SNES_LINESEARCH_CUBIC - cubic method
112178bcb3b5SPeter Brune 
112259405d5eSPeter Brune    Level: intermediate
112359405d5eSPeter Brune 
112459405d5eSPeter Brune .seealso: SNESLineSearchSetOrder()
112559405d5eSPeter Brune @*/
112659405d5eSPeter Brune 
1127b000cd8dSPeter Brune PetscErrorCode  SNESLineSearchGetOrder(SNESLineSearch linesearch,PetscInt *order)
112859405d5eSPeter Brune {
112959405d5eSPeter Brune   PetscFunctionBegin;
113059405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
113159405d5eSPeter Brune   PetscValidPointer(order, 2);
113259405d5eSPeter Brune   *order = linesearch->order;
113359405d5eSPeter Brune   PetscFunctionReturn(0);
113459405d5eSPeter Brune }
113559405d5eSPeter Brune 
113659405d5eSPeter Brune #undef __FUNCT__
113759405d5eSPeter Brune #define __FUNCT__ "SNESLineSearchSetOrder"
113859405d5eSPeter Brune /*@
113959405d5eSPeter Brune    SNESLineSearchSetOrder - Sets the line search damping paramter.
114059405d5eSPeter Brune 
114159405d5eSPeter Brune    Input Parameters:
114278bcb3b5SPeter Brune .  linesearch - linesearch context
114378bcb3b5SPeter Brune .  order - The damping parameter
114459405d5eSPeter Brune 
114559405d5eSPeter Brune    Level: intermediate
114659405d5eSPeter Brune 
114778bcb3b5SPeter Brune    Possible Values for order:
1148b000cd8dSPeter Brune +  SNES_LINESEARCH_ORDER_LINEAR - linear method
1149b000cd8dSPeter Brune .  SNES_LINESEARCH_ORDER_QUADRATIC - quadratic method
1150b000cd8dSPeter Brune -  SNES_LINESEARCH_ORDER_CUBIC - cubic method
115178bcb3b5SPeter Brune 
115259405d5eSPeter Brune    Notes:
115359405d5eSPeter Brune    Variable orders are supported by the following line searches:
115478bcb3b5SPeter Brune +  bt - cubic and quadratic
115578bcb3b5SPeter Brune -  cp - linear and quadratic
115659405d5eSPeter Brune 
115759405d5eSPeter Brune .seealso: SNESLineSearchGetOrder()
115859405d5eSPeter Brune @*/
1159b000cd8dSPeter Brune PetscErrorCode  SNESLineSearchSetOrder(SNESLineSearch linesearch,PetscInt order)
116059405d5eSPeter Brune {
116159405d5eSPeter Brune   PetscFunctionBegin;
116259405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
116359405d5eSPeter Brune   linesearch->order = order;
116459405d5eSPeter Brune   PetscFunctionReturn(0);
116559405d5eSPeter Brune }
116659405d5eSPeter Brune 
116759405d5eSPeter Brune #undef __FUNCT__
1168f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetNorms"
1169f40b411bSPeter Brune /*@
1170f1c6b773SPeter Brune    SNESLineSearchGetNorms - Gets the norms for for X, Y, and F.
1171f40b411bSPeter Brune 
1172f40b411bSPeter Brune    Input Parameters:
117378bcb3b5SPeter Brune .  linesearch - linesearch context
1174f40b411bSPeter Brune 
1175f40b411bSPeter Brune    Output Parameters:
1176f40b411bSPeter Brune +  xnorm - The norm of the current solution
1177f40b411bSPeter Brune .  fnorm - The norm of the current function
1178f40b411bSPeter Brune -  ynorm - The norm of the current update
1179f40b411bSPeter Brune 
1180cd7522eaSPeter Brune    Notes:
1181cd7522eaSPeter Brune    This function is mainly called from SNES implementations.
1182cd7522eaSPeter Brune 
118378bcb3b5SPeter Brune    Level: developer
1184f40b411bSPeter Brune 
1185f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms() SNESLineSearchGetVecs()
1186f40b411bSPeter Brune @*/
1187f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm)
1188bf7f4e0aSPeter Brune {
1189bf7f4e0aSPeter Brune   PetscFunctionBegin;
1190f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1191bf7f4e0aSPeter Brune   if (xnorm) {
1192bf7f4e0aSPeter Brune     *xnorm = linesearch->xnorm;
1193bf7f4e0aSPeter Brune   }
1194bf7f4e0aSPeter Brune   if (fnorm) {
1195bf7f4e0aSPeter Brune     *fnorm = linesearch->fnorm;
1196bf7f4e0aSPeter Brune   }
1197bf7f4e0aSPeter Brune   if (ynorm) {
1198bf7f4e0aSPeter Brune     *ynorm = linesearch->ynorm;
1199bf7f4e0aSPeter Brune   }
1200bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1201bf7f4e0aSPeter Brune }
1202bf7f4e0aSPeter Brune 
1203e7058c64SPeter Brune #undef __FUNCT__
1204f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetNorms"
1205f40b411bSPeter Brune /*@
1206f1c6b773SPeter Brune    SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F.
1207f40b411bSPeter Brune 
1208f40b411bSPeter Brune    Input Parameters:
120978bcb3b5SPeter Brune +  linesearch - linesearch context
1210f40b411bSPeter Brune .  xnorm - The norm of the current solution
1211f40b411bSPeter Brune .  fnorm - The norm of the current function
1212f40b411bSPeter Brune -  ynorm - The norm of the current update
1213f40b411bSPeter Brune 
121478bcb3b5SPeter Brune    Level: advanced
1215f40b411bSPeter Brune 
1216f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1217f40b411bSPeter Brune @*/
1218f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
12196a388c36SPeter Brune {
12206a388c36SPeter Brune   PetscFunctionBegin;
1221f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
12226a388c36SPeter Brune   linesearch->xnorm = xnorm;
12236a388c36SPeter Brune   linesearch->fnorm = fnorm;
12246a388c36SPeter Brune   linesearch->ynorm = ynorm;
12256a388c36SPeter Brune   PetscFunctionReturn(0);
12266a388c36SPeter Brune }
12276a388c36SPeter Brune 
12286a388c36SPeter Brune #undef __FUNCT__
1229f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchComputeNorms"
1230f40b411bSPeter Brune /*@
1231f1c6b773SPeter Brune    SNESLineSearchComputeNorms - Computes the norms of X, F, and Y.
1232f40b411bSPeter Brune 
1233f40b411bSPeter Brune    Input Parameters:
123478bcb3b5SPeter Brune .  linesearch - linesearch context
1235f40b411bSPeter Brune 
1236f40b411bSPeter Brune    Options Database Keys:
12378e557f58SPeter Brune .   -snes_linesearch_norms - turn norm computation on or off
1238f40b411bSPeter Brune 
1239f40b411bSPeter Brune    Level: intermediate
1240f40b411bSPeter Brune 
124178bcb3b5SPeter Brune .seealso: SNESLineSearchGetNorms, SNESLineSearchSetNorms(), SNESLineSearchSetComputeNorms()
1242f40b411bSPeter Brune @*/
1243f1c6b773SPeter Brune PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
12446a388c36SPeter Brune {
12456a388c36SPeter Brune   PetscErrorCode ierr;
12469bd66eb0SPeter Brune   SNES snes;
12476a388c36SPeter Brune   PetscFunctionBegin;
12486a388c36SPeter Brune   if (linesearch->norms) {
12499bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
1250f1c6b773SPeter Brune       ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
12519bd66eb0SPeter Brune       ierr = VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
12529bd66eb0SPeter Brune       ierr = VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
12539bd66eb0SPeter Brune       ierr = (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);CHKERRQ(ierr);
12549bd66eb0SPeter Brune     } else {
12556a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_func,   NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
12566a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_sol,    NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
12576a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
12586a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_func,     NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
12596a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_sol,      NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
12606a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_update,   NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
12616a388c36SPeter Brune     }
12629bd66eb0SPeter Brune   }
12636a388c36SPeter Brune   PetscFunctionReturn(0);
12646a388c36SPeter Brune }
12656a388c36SPeter Brune 
12666f263ca3SPeter Brune 
12676f263ca3SPeter Brune #undef __FUNCT__
12686f263ca3SPeter Brune #define __FUNCT__ "SNESLineSearchSetComputeNorms"
12696f263ca3SPeter Brune /*@
12706f263ca3SPeter Brune    SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search.
12716f263ca3SPeter Brune 
12726f263ca3SPeter Brune    Input Parameters:
127378bcb3b5SPeter Brune +  linesearch  - linesearch context
127478bcb3b5SPeter Brune -  flg  - indicates whether or not to compute norms
12756f263ca3SPeter Brune 
12766f263ca3SPeter Brune    Options Database Keys:
12778e557f58SPeter Brune .   -snes_linesearch_norms - turn norm computation on or off
12786f263ca3SPeter Brune 
12796f263ca3SPeter Brune    Notes:
12801a4f838cSPeter Brune    This is most relevant to the SNESLINESEARCHBASIC line search type.
12816f263ca3SPeter Brune 
12826f263ca3SPeter Brune    Level: intermediate
12836f263ca3SPeter Brune 
12841a4f838cSPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetNorms(), SNESLineSearchComputeNorms(), SNESLINESEARCHBASIC
12856f263ca3SPeter Brune @*/
12866f263ca3SPeter Brune PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg)
12876f263ca3SPeter Brune {
12886f263ca3SPeter Brune   PetscFunctionBegin;
12896f263ca3SPeter Brune   linesearch->norms = flg;
12906f263ca3SPeter Brune   PetscFunctionReturn(0);
12916f263ca3SPeter Brune }
12926f263ca3SPeter Brune 
12936a388c36SPeter Brune #undef __FUNCT__
1294f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetVecs"
1295f40b411bSPeter Brune /*@
1296f1c6b773SPeter Brune    SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context
1297f40b411bSPeter Brune 
1298f40b411bSPeter Brune    Input Parameters:
129978bcb3b5SPeter Brune .  linesearch - linesearch context
1300f40b411bSPeter Brune 
1301f40b411bSPeter Brune    Output Parameters:
1302f40b411bSPeter Brune +  X - The old solution
1303f40b411bSPeter Brune .  F - The old function
1304f40b411bSPeter Brune .  Y - The search direction
1305f40b411bSPeter Brune .  W - The new solution
1306f40b411bSPeter Brune -  G - The new function
1307f40b411bSPeter Brune 
130878bcb3b5SPeter Brune    Level: advanced
1309f40b411bSPeter Brune 
1310f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1311f40b411bSPeter Brune @*/
1312f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G) {
13136a388c36SPeter Brune   PetscFunctionBegin;
1314f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
13156a388c36SPeter Brune   if (X) {
13166a388c36SPeter Brune     PetscValidPointer(X, 2);
13176a388c36SPeter Brune     *X = linesearch->vec_sol;
13186a388c36SPeter Brune   }
13196a388c36SPeter Brune   if (F) {
13206a388c36SPeter Brune     PetscValidPointer(F, 3);
13216a388c36SPeter Brune     *F = linesearch->vec_func;
13226a388c36SPeter Brune   }
13236a388c36SPeter Brune   if (Y) {
13246a388c36SPeter Brune     PetscValidPointer(Y, 4);
13256a388c36SPeter Brune     *Y = linesearch->vec_update;
13266a388c36SPeter Brune   }
13276a388c36SPeter Brune   if (W) {
13286a388c36SPeter Brune     PetscValidPointer(W, 5);
13296a388c36SPeter Brune     *W = linesearch->vec_sol_new;
13306a388c36SPeter Brune   }
13316a388c36SPeter Brune   if (G) {
13326a388c36SPeter Brune     PetscValidPointer(G, 6);
13336a388c36SPeter Brune     *G = linesearch->vec_func_new;
13346a388c36SPeter Brune   }
13356a388c36SPeter Brune 
13366a388c36SPeter Brune   PetscFunctionReturn(0);
13376a388c36SPeter Brune }
13386a388c36SPeter Brune 
13396a388c36SPeter Brune #undef __FUNCT__
1340f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetVecs"
1341f40b411bSPeter Brune /*@
1342f1c6b773SPeter Brune    SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context
1343f40b411bSPeter Brune 
1344f40b411bSPeter Brune    Input Parameters:
134578bcb3b5SPeter Brune +  linesearch - linesearch context
1346f40b411bSPeter Brune .  X - The old solution
1347f40b411bSPeter Brune .  F - The old function
1348f40b411bSPeter Brune .  Y - The search direction
1349f40b411bSPeter Brune .  W - The new solution
1350f40b411bSPeter Brune -  G - The new function
1351f40b411bSPeter Brune 
135278bcb3b5SPeter Brune    Level: advanced
1353f40b411bSPeter Brune 
1354f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms(), SNESLineSearchGetVecs()
1355f40b411bSPeter Brune @*/
1356f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G) {
13576a388c36SPeter Brune   PetscFunctionBegin;
1358f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
13596a388c36SPeter Brune   if (X) {
13606a388c36SPeter Brune     PetscValidHeaderSpecific(X,VEC_CLASSID,2);
13616a388c36SPeter Brune     linesearch->vec_sol = X;
13626a388c36SPeter Brune   }
13636a388c36SPeter Brune   if (F) {
13646a388c36SPeter Brune     PetscValidHeaderSpecific(F,VEC_CLASSID,3);
13656a388c36SPeter Brune     linesearch->vec_func = F;
13666a388c36SPeter Brune   }
13676a388c36SPeter Brune   if (Y) {
13686a388c36SPeter Brune     PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
13696a388c36SPeter Brune     linesearch->vec_update = Y;
13706a388c36SPeter Brune   }
13716a388c36SPeter Brune   if (W) {
13726a388c36SPeter Brune     PetscValidHeaderSpecific(W,VEC_CLASSID,5);
13736a388c36SPeter Brune     linesearch->vec_sol_new = W;
13746a388c36SPeter Brune   }
13756a388c36SPeter Brune   if (G) {
13766a388c36SPeter Brune     PetscValidHeaderSpecific(G,VEC_CLASSID,6);
13776a388c36SPeter Brune     linesearch->vec_func_new = G;
13786a388c36SPeter Brune   }
13796a388c36SPeter Brune 
13806a388c36SPeter Brune   PetscFunctionReturn(0);
13816a388c36SPeter Brune }
13826a388c36SPeter Brune 
13836a388c36SPeter Brune #undef __FUNCT__
1384f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchAppendOptionsPrefix"
1385e7058c64SPeter Brune /*@C
1386f1c6b773SPeter Brune    SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1387e7058c64SPeter Brune    SNES options in the database.
1388e7058c64SPeter Brune 
1389cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
1390e7058c64SPeter Brune 
1391e7058c64SPeter Brune    Input Parameters:
1392e7058c64SPeter Brune +  snes - the SNES context
1393e7058c64SPeter Brune -  prefix - the prefix to prepend to all option names
1394e7058c64SPeter Brune 
1395e7058c64SPeter Brune    Notes:
1396e7058c64SPeter Brune    A hyphen (-) must NOT be given at the beginning of the prefix name.
1397e7058c64SPeter Brune    The first character of all runtime options is AUTOMATICALLY the hyphen.
1398e7058c64SPeter Brune 
1399e7058c64SPeter Brune    Level: advanced
1400e7058c64SPeter Brune 
1401f1c6b773SPeter Brune .keywords: SNESLineSearch, append, options, prefix, database
1402e7058c64SPeter Brune 
1403e7058c64SPeter Brune .seealso: SNESGetOptionsPrefix()
1404e7058c64SPeter Brune @*/
1405f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[])
1406e7058c64SPeter Brune {
1407e7058c64SPeter Brune   PetscErrorCode ierr;
1408e7058c64SPeter Brune 
1409e7058c64SPeter Brune   PetscFunctionBegin;
1410f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1411e7058c64SPeter Brune   ierr = PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1412e7058c64SPeter Brune   PetscFunctionReturn(0);
1413e7058c64SPeter Brune }
1414e7058c64SPeter Brune 
1415e7058c64SPeter Brune #undef __FUNCT__
1416f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetOptionsPrefix"
1417e7058c64SPeter Brune /*@C
1418f1c6b773SPeter Brune    SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all
1419f1c6b773SPeter Brune    SNESLineSearch options in the database.
1420e7058c64SPeter Brune 
1421e7058c64SPeter Brune    Not Collective
1422e7058c64SPeter Brune 
1423e7058c64SPeter Brune    Input Parameter:
1424cd7522eaSPeter Brune .  linesearch - the SNESLineSearch context
1425e7058c64SPeter Brune 
1426e7058c64SPeter Brune    Output Parameter:
1427e7058c64SPeter Brune .  prefix - pointer to the prefix string used
1428e7058c64SPeter Brune 
14298e557f58SPeter Brune    Notes:
14308e557f58SPeter Brune    On the fortran side, the user should pass in a string 'prefix' of
1431e7058c64SPeter Brune    sufficient length to hold the prefix.
1432e7058c64SPeter Brune 
1433e7058c64SPeter Brune    Level: advanced
1434e7058c64SPeter Brune 
1435f1c6b773SPeter Brune .keywords: SNESLineSearch, get, options, prefix, database
1436e7058c64SPeter Brune 
1437e7058c64SPeter Brune .seealso: SNESAppendOptionsPrefix()
1438e7058c64SPeter Brune @*/
1439f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[])
1440e7058c64SPeter Brune {
1441e7058c64SPeter Brune   PetscErrorCode ierr;
1442e7058c64SPeter Brune 
1443e7058c64SPeter Brune   PetscFunctionBegin;
1444f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1445e7058c64SPeter Brune   ierr = PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1446e7058c64SPeter Brune   PetscFunctionReturn(0);
1447e7058c64SPeter Brune }
1448bf7f4e0aSPeter Brune 
1449bf7f4e0aSPeter Brune #undef __FUNCT__
1450f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetWork"
1451f40b411bSPeter Brune /*@
1452f1c6b773SPeter Brune    SNESLineSearchGetWork - Gets work vectors for the line search.
1453f40b411bSPeter Brune 
1454f40b411bSPeter Brune    Input Parameter:
1455f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
1456f40b411bSPeter Brune -  nwork - the number of work vectors
1457f40b411bSPeter Brune 
1458f40b411bSPeter Brune    Level: developer
1459f40b411bSPeter Brune 
1460cd7522eaSPeter Brune    Notes:
1461cd7522eaSPeter Brune    This is typically called at the beginning of a SNESLineSearch or SNESLineSearchShell implementation.
1462cd7522eaSPeter Brune 
1463f1c6b773SPeter Brune .keywords: SNESLineSearch, work, vector
1464f40b411bSPeter Brune 
1465f40b411bSPeter Brune .seealso: SNESDefaultGetWork()
1466f40b411bSPeter Brune @*/
1467f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetWork(SNESLineSearch linesearch, PetscInt nwork)
1468bf7f4e0aSPeter Brune {
1469bf7f4e0aSPeter Brune   PetscErrorCode ierr;
1470bf7f4e0aSPeter Brune   PetscFunctionBegin;
1471bf7f4e0aSPeter Brune   if (linesearch->vec_sol) {
1472bf7f4e0aSPeter Brune     ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr);
1473bf7f4e0aSPeter Brune   } else {
1474bf7f4e0aSPeter Brune     SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
1475bf7f4e0aSPeter Brune   }
1476bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1477bf7f4e0aSPeter Brune }
1478bf7f4e0aSPeter Brune 
14796a388c36SPeter Brune 
1480bf7f4e0aSPeter Brune #undef __FUNCT__
1481f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetSuccess"
1482f40b411bSPeter Brune /*@
1483f1c6b773SPeter Brune    SNESLineSearchGetSuccess - Gets the success/failure status of the last line search application
1484f40b411bSPeter Brune 
1485f40b411bSPeter Brune    Input Parameters:
148678bcb3b5SPeter Brune .  linesearch - linesearch context
1487f40b411bSPeter Brune 
1488f40b411bSPeter Brune    Output Parameters:
14898e557f58SPeter Brune .  success - The success or failure status
1490f40b411bSPeter Brune 
1491cd7522eaSPeter Brune    Notes:
1492cd7522eaSPeter Brune    This is typically called after SNESLineSearchApply in order to determine if the line-search failed
1493cd7522eaSPeter Brune    (and set the SNES convergence accordingly).
1494cd7522eaSPeter Brune 
1495f40b411bSPeter Brune    Level: intermediate
1496f40b411bSPeter Brune 
1497f1c6b773SPeter Brune .seealso: SNESLineSearchSetSuccess()
1498f40b411bSPeter Brune @*/
1499f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetSuccess(SNESLineSearch linesearch, PetscBool *success)
1500bf7f4e0aSPeter Brune {
1501bf7f4e0aSPeter Brune   PetscFunctionBegin;
1502f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
15036a388c36SPeter Brune   PetscValidPointer(success, 2);
1504bf7f4e0aSPeter Brune   if (success) {
1505bf7f4e0aSPeter Brune     *success = linesearch->success;
1506bf7f4e0aSPeter Brune   }
1507bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1508bf7f4e0aSPeter Brune }
1509bf7f4e0aSPeter Brune 
1510bf7f4e0aSPeter Brune #undef __FUNCT__
1511f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetSuccess"
1512f40b411bSPeter Brune /*@
1513f1c6b773SPeter Brune    SNESLineSearchSetSuccess - Sets the success/failure status of the last line search application
1514f40b411bSPeter Brune 
1515f40b411bSPeter Brune    Input Parameters:
151678bcb3b5SPeter Brune +  linesearch - linesearch context
15178e557f58SPeter Brune -  success - The success or failure status
1518f40b411bSPeter Brune 
1519cd7522eaSPeter Brune    Notes:
1520cd7522eaSPeter Brune    This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set
1521cd7522eaSPeter Brune    the success or failure of the line search method.
1522cd7522eaSPeter Brune 
152378bcb3b5SPeter Brune    Level: developer
1524f40b411bSPeter Brune 
1525f1c6b773SPeter Brune .seealso: SNESLineSearchGetSuccess()
1526f40b411bSPeter Brune @*/
1527f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetSuccess(SNESLineSearch linesearch, PetscBool success)
15286a388c36SPeter Brune {
1529f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
15306a388c36SPeter Brune   PetscFunctionBegin;
15316a388c36SPeter Brune   linesearch->success = success;
15326a388c36SPeter Brune   PetscFunctionReturn(0);
15336a388c36SPeter Brune }
15346a388c36SPeter Brune 
15356a388c36SPeter Brune #undef __FUNCT__
1536f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetVIFunctions"
15379bd66eb0SPeter Brune /*@C
1538f1c6b773SPeter Brune    SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.
15399bd66eb0SPeter Brune 
15409bd66eb0SPeter Brune    Input Parameters:
15419bd66eb0SPeter Brune +  snes - nonlinear context obtained from SNESCreate()
15429bd66eb0SPeter Brune .  projectfunc - function for projecting the function to the bounds
15439bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
15449bd66eb0SPeter Brune 
15459bd66eb0SPeter Brune    Logically Collective on SNES
15469bd66eb0SPeter Brune 
15479bd66eb0SPeter Brune    Calling sequence of projectfunc:
15489bd66eb0SPeter Brune .vb
15499bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X)
15509bd66eb0SPeter Brune .ve
15519bd66eb0SPeter Brune 
15529bd66eb0SPeter Brune     Input parameters for projectfunc:
15539bd66eb0SPeter Brune +   snes - nonlinear context
15549bd66eb0SPeter Brune -   X - current solution
15559bd66eb0SPeter Brune 
1556cd7522eaSPeter Brune     Output parameters for projectfunc:
15579bd66eb0SPeter Brune .   X - Projected solution
15589bd66eb0SPeter Brune 
15599bd66eb0SPeter Brune    Calling sequence of normfunc:
15609bd66eb0SPeter Brune .vb
15619bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm)
15629bd66eb0SPeter Brune .ve
15639bd66eb0SPeter Brune 
1564cd7522eaSPeter Brune     Input parameters for normfunc:
15659bd66eb0SPeter Brune +   snes - nonlinear context
15669bd66eb0SPeter Brune .   X - current solution
15679bd66eb0SPeter Brune -   F - current residual
15689bd66eb0SPeter Brune 
1569cd7522eaSPeter Brune     Output parameters for normfunc:
15709bd66eb0SPeter Brune .   fnorm - VI-specific norm of the function
15719bd66eb0SPeter Brune 
1572cd7522eaSPeter Brune     Notes:
1573cd7522eaSPeter Brune     The VI solvers require projection of the solution to the feasible set.  projectfunc should implement this.
1574cd7522eaSPeter Brune 
1575cd7522eaSPeter Brune     The VI solvers require special evaluation of the function norm such that the norm is only calculated
1576cd7522eaSPeter Brune     on the inactive set.  This should be implemented by normfunc.
15779bd66eb0SPeter Brune 
15789bd66eb0SPeter Brune     Level: developer
15799bd66eb0SPeter Brune 
15809bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, set, line search
15819bd66eb0SPeter Brune 
1582f1c6b773SPeter Brune .seealso: SNESLineSearchGetVIFunctions(), SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck()
15839bd66eb0SPeter Brune @*/
1584f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc)
15859bd66eb0SPeter Brune {
15869bd66eb0SPeter Brune   PetscFunctionBegin;
1587f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
15889bd66eb0SPeter Brune   if (projectfunc) linesearch->ops->viproject = projectfunc;
15899bd66eb0SPeter Brune   if (normfunc) linesearch->ops->vinorm = normfunc;
15909bd66eb0SPeter Brune   PetscFunctionReturn(0);
15919bd66eb0SPeter Brune }
15929bd66eb0SPeter Brune 
15939bd66eb0SPeter Brune #undef __FUNCT__
1594f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetVIFunctions"
15959bd66eb0SPeter Brune /*@C
1596f1c6b773SPeter Brune    SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.
15979bd66eb0SPeter Brune 
15989bd66eb0SPeter Brune    Input Parameters:
15999bd66eb0SPeter Brune .  snes - nonlinear context obtained from SNESCreate()
16009bd66eb0SPeter Brune 
16019bd66eb0SPeter Brune    Output Parameters:
16029bd66eb0SPeter Brune +  projectfunc - function for projecting the function to the bounds
16039bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
16049bd66eb0SPeter Brune 
16059bd66eb0SPeter Brune    Logically Collective on SNES
16069bd66eb0SPeter Brune 
16079bd66eb0SPeter Brune     Level: developer
16089bd66eb0SPeter Brune 
16099bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, get, line search
16109bd66eb0SPeter Brune 
1611f1c6b773SPeter Brune .seealso: SNESLineSearchSetVIFunctions(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
16129bd66eb0SPeter Brune @*/
1613f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc)
16149bd66eb0SPeter Brune {
16159bd66eb0SPeter Brune   PetscFunctionBegin;
16169bd66eb0SPeter Brune   if (projectfunc) *projectfunc = linesearch->ops->viproject;
16179bd66eb0SPeter Brune   if (normfunc) *normfunc = linesearch->ops->vinorm;
16189bd66eb0SPeter Brune   PetscFunctionReturn(0);
16199bd66eb0SPeter Brune }
16209bd66eb0SPeter Brune 
16219bd66eb0SPeter Brune #undef __FUNCT__
1622f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchRegister"
1623bf7f4e0aSPeter Brune /*@C
1624f1c6b773SPeter Brune   SNESLineSearchRegister - See SNESLineSearchRegisterDynamic()
1625bf7f4e0aSPeter Brune 
1626bf7f4e0aSPeter Brune   Level: advanced
1627bf7f4e0aSPeter Brune @*/
1628f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNESLineSearch))
1629bf7f4e0aSPeter Brune {
1630bf7f4e0aSPeter Brune   char           fullname[PETSC_MAX_PATH_LEN];
1631bf7f4e0aSPeter Brune   PetscErrorCode ierr;
1632bf7f4e0aSPeter Brune 
1633bf7f4e0aSPeter Brune   PetscFunctionBegin;
1634bf7f4e0aSPeter Brune   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
1635f1c6b773SPeter Brune   ierr = PetscFListAdd(&SNESLineSearchList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
1636bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1637bf7f4e0aSPeter Brune }
1638