xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision 6f263ca3b312f6197d13b46ac2297b18ba3ebe1a)
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:
20f40b411bSPeter Brune .  outlinesearch - the line search instance.
21f40b411bSPeter Brune 
22f40b411bSPeter 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 /*@
71f1c6b773SPeter Brune    SNESLineSearchSetUp - Prepares the line search for being applied.
72f40b411bSPeter Brune 
73cd7522eaSPeter Brune    Collective on SNESLineSearch
74f40b411bSPeter Brune 
75f40b411bSPeter Brune    Input Parameters:
76f40b411bSPeter Brune .  linesearch - The LineSearch instance.
77f40b411bSPeter Brune 
78cd7522eaSPeter Brune    Notes:
79cd7522eaSPeter Brune    For most cases, this needn't be called outside of SNESLineSearchApply().
80cd7522eaSPeter Brune    The only current case where this is called outside of this is for the VI
81cd7522eaSPeter Brune    solvers, which modifies the solution and work vectors before the first call
82cd7522eaSPeter Brune    of SNESLineSearchApply, requiring the SNESLineSearch work vectors to be
83cd7522eaSPeter Brune    allocated upfront.
84cd7522eaSPeter Brune 
85cd7522eaSPeter Brune 
86f40b411bSPeter Brune    Level: Intermediate
87f40b411bSPeter Brune 
88f1c6b773SPeter Brune    .keywords: SNESLineSearch, SetUp
89f40b411bSPeter Brune 
90f1c6b773SPeter Brune    .seealso: SNESLineSearchReset()
91f40b411bSPeter Brune @*/
92f40b411bSPeter Brune 
93f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch) {
94bf7f4e0aSPeter Brune   PetscErrorCode ierr;
95bf7f4e0aSPeter Brune   PetscFunctionBegin;
96bf7f4e0aSPeter Brune   if (!((PetscObject)linesearch)->type_name) {
97f1c6b773SPeter Brune     ierr = SNESLineSearchSetType(linesearch,SNES_LINESEARCH_BASIC);CHKERRQ(ierr);
98bf7f4e0aSPeter Brune   }
99bf7f4e0aSPeter Brune   if (!linesearch->setupcalled) {
1009bd66eb0SPeter Brune     if (!linesearch->vec_sol_new) {
101bf7f4e0aSPeter Brune       ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new);CHKERRQ(ierr);
1029bd66eb0SPeter Brune     }
1039bd66eb0SPeter Brune     if (!linesearch->vec_func_new) {
1049bd66eb0SPeter Brune       ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new);CHKERRQ(ierr);
1059bd66eb0SPeter Brune     }
106bf7f4e0aSPeter Brune     if (linesearch->ops->setup) {
107bf7f4e0aSPeter Brune       ierr = (*linesearch->ops->setup)(linesearch);CHKERRQ(ierr);
108bf7f4e0aSPeter Brune     }
109bf7f4e0aSPeter Brune     linesearch->lambda = linesearch->damping;
110bf7f4e0aSPeter Brune     linesearch->setupcalled = PETSC_TRUE;
111bf7f4e0aSPeter Brune   }
112bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
113bf7f4e0aSPeter Brune }
114bf7f4e0aSPeter Brune 
115bf7f4e0aSPeter Brune #undef __FUNCT__
116f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchReset"
117f40b411bSPeter Brune 
118f40b411bSPeter Brune /*@
119cd7522eaSPeter Brune    SNESLineSearchReset - Tears down the structures required for application of the SNESLineSearch
120f40b411bSPeter Brune 
121f1c6b773SPeter Brune    Collective on SNESLineSearch
122f40b411bSPeter Brune 
123f40b411bSPeter Brune    Input Parameters:
124f40b411bSPeter Brune .  linesearch - The LineSearch instance.
125f40b411bSPeter Brune 
126f40b411bSPeter Brune    Level: Intermediate
127f40b411bSPeter Brune 
128cd7522eaSPeter Brune    .keywords: SNESLineSearch, Reset
129f40b411bSPeter Brune 
130f1c6b773SPeter Brune    .seealso: SNESLineSearchSetUp()
131f40b411bSPeter Brune @*/
132f40b411bSPeter Brune 
133f1c6b773SPeter Brune PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch) {
134bf7f4e0aSPeter Brune   PetscErrorCode ierr;
135bf7f4e0aSPeter Brune   PetscFunctionBegin;
136bf7f4e0aSPeter Brune   if (linesearch->ops->reset) {
137bf7f4e0aSPeter Brune     (*linesearch->ops->reset)(linesearch);
138bf7f4e0aSPeter Brune   }
139bf7f4e0aSPeter Brune   ierr = VecDestroy(&linesearch->vec_sol_new);CHKERRQ(ierr);
140bf7f4e0aSPeter Brune   ierr = VecDestroy(&linesearch->vec_func_new);CHKERRQ(ierr);
141bf7f4e0aSPeter Brune 
142bf7f4e0aSPeter Brune   ierr = VecDestroyVecs(linesearch->nwork, &linesearch->work);CHKERRQ(ierr);
143bf7f4e0aSPeter Brune   linesearch->nwork = 0;
144bf7f4e0aSPeter Brune   linesearch->setupcalled = PETSC_FALSE;
145bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
146bf7f4e0aSPeter Brune }
147bf7f4e0aSPeter Brune 
14886d74e61SPeter Brune 
14986d74e61SPeter Brune #undef __FUNCT__
150f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetPreCheck"
15186d74e61SPeter Brune /*@C
152f1c6b773SPeter Brune    SNESLineSearchSetPreCheck - Sets a pre-check function for the line search routine.
15386d74e61SPeter Brune 
154f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
15586d74e61SPeter Brune 
15686d74e61SPeter Brune    Input Parameters:
157f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
15886d74e61SPeter Brune .  func       - [optional] function evaluation routine
15986d74e61SPeter Brune -  ctx        - [optional] user-defined context for private data for the
16086d74e61SPeter Brune                 function evaluation routine (may be PETSC_NULL)
16186d74e61SPeter Brune 
16286d74e61SPeter Brune    Calling sequence of func:
163f1c6b773SPeter Brune $    func (SNESLineSearch snes,Vec x,Vec y, PetscBool *changed);
16486d74e61SPeter Brune 
16586d74e61SPeter Brune +  x - solution vector
16686d74e61SPeter Brune .  y - search direction vector
167cd7522eaSPeter Brune -  changed - flag to indicate the precheck changed x or y.
16886d74e61SPeter Brune 
16986d74e61SPeter Brune    Level: intermediate
17086d74e61SPeter Brune 
17186d74e61SPeter Brune .keywords: set, linesearch, pre-check
17286d74e61SPeter Brune 
173f1c6b773SPeter Brune .seealso: SNESLineSearchSetPostCheck()
17486d74e61SPeter Brune @*/
175f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetPreCheck(SNESLineSearch linesearch, SNESLineSearchPreCheckFunc func,void *ctx)
17686d74e61SPeter Brune {
1779bd66eb0SPeter Brune   PetscFunctionBegin;
178f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
17986d74e61SPeter Brune   if (func) linesearch->ops->precheckstep = func;
18086d74e61SPeter Brune   if (ctx) linesearch->precheckctx = ctx;
18186d74e61SPeter Brune   PetscFunctionReturn(0);
18286d74e61SPeter Brune }
18386d74e61SPeter Brune 
18486d74e61SPeter Brune 
18586d74e61SPeter Brune #undef __FUNCT__
186f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetPreCheck"
18786d74e61SPeter Brune /*@C
188cd7522eaSPeter Brune    SNESLineSearchSetPreCheck - Gets the pre-check function for the line search routine.
18986d74e61SPeter Brune 
19086d74e61SPeter Brune    Input Parameters:
191f1c6b773SPeter Brune .  linesearch - the SNESLineSearch context
19286d74e61SPeter Brune 
19386d74e61SPeter Brune    Output Parameters:
19486d74e61SPeter Brune +  func       - [optional] function evaluation routine
19586d74e61SPeter Brune -  ctx        - [optional] user-defined context for private data for the
19686d74e61SPeter Brune                 function evaluation routine (may be PETSC_NULL)
19786d74e61SPeter Brune 
19886d74e61SPeter Brune    Level: intermediate
19986d74e61SPeter Brune 
20086d74e61SPeter Brune .keywords: get, linesearch, pre-check
20186d74e61SPeter Brune 
202f1c6b773SPeter Brune .seealso: SNESLineSearchGetPostCheck(), SNESLineSearchSetPreCheck()
20386d74e61SPeter Brune @*/
204f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetPreCheck(SNESLineSearch linesearch, SNESLineSearchPreCheckFunc *func,void **ctx)
20586d74e61SPeter Brune {
2069bd66eb0SPeter Brune   PetscFunctionBegin;
207f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
20886d74e61SPeter Brune   if (func) *func = linesearch->ops->precheckstep;
20986d74e61SPeter Brune   if (ctx) *ctx = linesearch->precheckctx;
21086d74e61SPeter Brune   PetscFunctionReturn(0);
21186d74e61SPeter Brune }
21286d74e61SPeter Brune 
21386d74e61SPeter Brune 
21486d74e61SPeter Brune #undef __FUNCT__
215f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetPostCheck"
21686d74e61SPeter Brune /*@C
217f1c6b773SPeter Brune    SNESLineSearchSetPostCheck - Sets a post-check function for the line search routine.
21886d74e61SPeter Brune 
219f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
22086d74e61SPeter Brune 
22186d74e61SPeter Brune    Input Parameters:
222f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
22386d74e61SPeter Brune .  func       - [optional] function evaluation routine
22486d74e61SPeter Brune -  ctx        - [optional] user-defined context for private data for the
22586d74e61SPeter Brune                 function evaluation routine (may be PETSC_NULL)
22686d74e61SPeter Brune 
22786d74e61SPeter Brune    Calling sequence of func:
228f1c6b773SPeter Brune $    func (SNESLineSearch linesearch,Vec x, Vec w, Vec y, PetscBool *changed_w, *changed_y);
22986d74e61SPeter Brune 
23086d74e61SPeter Brune +  x - old solution vector
23186d74e61SPeter Brune .  y - search direction vector
23286d74e61SPeter Brune .  w - new solution vector
23386d74e61SPeter Brune .  changed_y - indicates that the line search changed y.
23486d74e61SPeter Brune .  changed_w - indicates that the line search changed w.
23586d74e61SPeter Brune 
23686d74e61SPeter Brune    Level: intermediate
23786d74e61SPeter Brune 
23886d74e61SPeter Brune .keywords: set, linesearch, post-check
23986d74e61SPeter Brune 
240f1c6b773SPeter Brune .seealso: SNESLineSearchSetPreCheck()
24186d74e61SPeter Brune @*/
242f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetPostCheck(SNESLineSearch linesearch, SNESLineSearchPostCheckFunc func,void *ctx)
24386d74e61SPeter Brune {
24486d74e61SPeter Brune   PetscFunctionBegin;
245f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
24686d74e61SPeter Brune   if (func) linesearch->ops->postcheckstep = func;
24786d74e61SPeter Brune   if (ctx) linesearch->postcheckctx = ctx;
24886d74e61SPeter Brune   PetscFunctionReturn(0);
24986d74e61SPeter Brune }
25086d74e61SPeter Brune 
25186d74e61SPeter Brune 
25286d74e61SPeter Brune #undef __FUNCT__
253f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetPostCheck"
25486d74e61SPeter Brune /*@C
255f1c6b773SPeter Brune    SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine.
25686d74e61SPeter Brune 
25786d74e61SPeter Brune    Input Parameters:
258f1c6b773SPeter Brune .  linesearch - the SNESLineSearch context
25986d74e61SPeter Brune 
26086d74e61SPeter Brune    Output Parameters:
26186d74e61SPeter Brune +  func       - [optional] function evaluation routine
26286d74e61SPeter Brune -  ctx        - [optional] user-defined context for private data for the
26386d74e61SPeter Brune                 function evaluation routine (may be PETSC_NULL)
26486d74e61SPeter Brune 
26586d74e61SPeter Brune    Level: intermediate
26686d74e61SPeter Brune 
26786d74e61SPeter Brune .keywords: get, linesearch, post-check
26886d74e61SPeter Brune 
269f1c6b773SPeter Brune .seealso: SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck()
27086d74e61SPeter Brune @*/
271f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetPostCheck(SNESLineSearch linesearch, SNESLineSearchPostCheckFunc *func,void **ctx)
27286d74e61SPeter Brune {
2739bd66eb0SPeter Brune   PetscFunctionBegin;
274f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
27586d74e61SPeter Brune   if (func) *func = linesearch->ops->postcheckstep;
27686d74e61SPeter Brune   if (ctx) *ctx = linesearch->postcheckctx;
27786d74e61SPeter Brune   PetscFunctionReturn(0);
27886d74e61SPeter Brune }
27986d74e61SPeter Brune 
28086d74e61SPeter Brune 
281bf7f4e0aSPeter Brune #undef __FUNCT__
282f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPreCheck"
283f40b411bSPeter Brune /*@
284f1c6b773SPeter Brune    SNESLineSearchPreCheck - Prepares the line search for being applied.
285f40b411bSPeter Brune 
286cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
287f40b411bSPeter Brune 
288f40b411bSPeter Brune    Input Parameters:
289f40b411bSPeter Brune .  linesearch - The linesearch instance.
290f40b411bSPeter Brune 
291f40b411bSPeter Brune    Output Parameters:
292f40b411bSPeter Brune .  changed - Indicator if the pre-check has changed anything.
293f40b411bSPeter Brune 
294f40b411bSPeter Brune    Level: Beginner
295f40b411bSPeter Brune 
296f1c6b773SPeter Brune    .keywords: SNESLineSearch, Create
297f40b411bSPeter Brune 
298f1c6b773SPeter Brune    .seealso: SNESLineSearchPostCheck()
299f40b411bSPeter Brune @*/
300f1c6b773SPeter Brune PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch, PetscBool * changed)
301bf7f4e0aSPeter Brune {
302bf7f4e0aSPeter Brune   PetscErrorCode ierr;
303bf7f4e0aSPeter Brune   PetscFunctionBegin;
304bf7f4e0aSPeter Brune   *changed = PETSC_FALSE;
305bf7f4e0aSPeter Brune   if (linesearch->ops->precheckstep) {
306c87759e9SPeter Brune     ierr = (*linesearch->ops->precheckstep)(linesearch, linesearch->vec_sol, linesearch->vec_update, changed, linesearch->precheckctx);CHKERRQ(ierr);
307bf7f4e0aSPeter Brune   }
308bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
309bf7f4e0aSPeter Brune }
310bf7f4e0aSPeter Brune 
311bf7f4e0aSPeter Brune #undef __FUNCT__
312f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPostCheck"
313f40b411bSPeter Brune /*@
314f1c6b773SPeter Brune    SNESLineSearchPostCheck - Prepares the line search for being applied.
315f40b411bSPeter Brune 
316cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
317f40b411bSPeter Brune 
318f40b411bSPeter Brune    Input Parameters:
319f40b411bSPeter Brune .  linesearch - The linesearch instance.
320f40b411bSPeter Brune 
321f40b411bSPeter Brune    Output Parameters:
32268fb48e8SJed Brown +  changed_Y - Indicator if the direction has been changed.
32368fb48e8SJed Brown -  changed_W - Indicator if the solution has been changed.
324f40b411bSPeter Brune 
325f40b411bSPeter Brune    Level: Intermediate
326f40b411bSPeter Brune 
327f1c6b773SPeter Brune    .keywords: SNESLineSearch, Create
328f40b411bSPeter Brune 
329f1c6b773SPeter Brune    .seealso: SNESLineSearchPreCheck()
330f40b411bSPeter Brune @*/
331f1c6b773SPeter Brune PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch, PetscBool * changed_Y, PetscBool * changed_W)
332bf7f4e0aSPeter Brune {
333bf7f4e0aSPeter Brune   PetscErrorCode ierr;
334bf7f4e0aSPeter Brune   PetscFunctionBegin;
335bf7f4e0aSPeter Brune   *changed_Y = PETSC_FALSE;
336bf7f4e0aSPeter Brune   *changed_W = PETSC_FALSE;
337bf7f4e0aSPeter Brune   if (linesearch->ops->postcheckstep) {
338c87759e9SPeter Brune     ierr = (*linesearch->ops->postcheckstep)(linesearch, linesearch->vec_sol, linesearch->vec_update, linesearch->vec_sol_new, changed_Y, changed_W, linesearch->postcheckctx);CHKERRQ(ierr);
33986d74e61SPeter Brune   }
34086d74e61SPeter Brune   PetscFunctionReturn(0);
34186d74e61SPeter Brune }
34286d74e61SPeter Brune 
34386d74e61SPeter Brune 
34486d74e61SPeter Brune #undef __FUNCT__
345f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPreCheckPicard"
34686d74e61SPeter Brune /*@C
34786d74e61SPeter Brune    SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration
34886d74e61SPeter Brune 
349cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
35086d74e61SPeter Brune 
35186d74e61SPeter Brune    Input Arguments:
35286d74e61SPeter Brune +  linesearch - linesearch context
35386d74e61SPeter Brune .  X - base state for this step
35486d74e61SPeter Brune .  Y - initial correction
35586d74e61SPeter Brune 
35686d74e61SPeter Brune    Output Arguments:
35786d74e61SPeter Brune +  Y - correction, possibly modified
35886d74e61SPeter Brune -  changed - flag indicating that Y was modified
35986d74e61SPeter Brune 
36086d74e61SPeter Brune    Options Database Key:
361cd7522eaSPeter Brune +  -snes_linesearch_precheck_picard - activate this routine
362cd7522eaSPeter Brune -  -snes_linesearch_precheck_picard_angle - angle
36386d74e61SPeter Brune 
36486d74e61SPeter Brune    Level: advanced
36586d74e61SPeter Brune 
36686d74e61SPeter Brune    Notes:
36786d74e61SPeter Brune    This function should be passed to SNESLineSearchSetPreCheck()
36886d74e61SPeter Brune 
36986d74e61SPeter Brune    The justification for this method involves the linear convergence of a Picard iteration
37086d74e61SPeter Brune    so the Picard linearization should be provided in place of the "Jacobian". This correction
37186d74e61SPeter Brune    is generally not useful when using a Newton linearization.
37286d74e61SPeter Brune 
37386d74e61SPeter Brune    Reference:
37486d74e61SPeter Brune    Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology.
37586d74e61SPeter Brune 
37686d74e61SPeter Brune .seealso: SNESLineSearchSetPreCheck()
37786d74e61SPeter Brune @*/
378f1c6b773SPeter Brune PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed,void *ctx)
37986d74e61SPeter Brune {
38086d74e61SPeter Brune   PetscErrorCode ierr;
38186d74e61SPeter Brune   PetscReal      angle = *(PetscReal*)linesearch->precheckctx;
38286d74e61SPeter Brune   Vec            Ylast;
38386d74e61SPeter Brune   PetscScalar    dot;
384c87759e9SPeter Brune 
38586d74e61SPeter Brune   PetscInt       iter;
38686d74e61SPeter Brune   PetscReal      ynorm,ylastnorm,theta,angle_radians;
38786d74e61SPeter Brune   SNES           snes;
38886d74e61SPeter Brune 
38986d74e61SPeter Brune   PetscFunctionBegin;
390f1c6b773SPeter Brune   ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
39186d74e61SPeter Brune   ierr = PetscObjectQuery((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject*)&Ylast);CHKERRQ(ierr);
39286d74e61SPeter Brune   if (!Ylast) {
39386d74e61SPeter Brune     ierr = VecDuplicate(Y,&Ylast);CHKERRQ(ierr);
39486d74e61SPeter Brune     ierr = PetscObjectCompose((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject)Ylast);CHKERRQ(ierr);
39586d74e61SPeter Brune     ierr = PetscObjectDereference((PetscObject)Ylast);CHKERRQ(ierr);
39686d74e61SPeter Brune   }
39786d74e61SPeter Brune   ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr);
39886d74e61SPeter Brune   if (iter < 2) {
39986d74e61SPeter Brune     ierr = VecCopy(Y,Ylast);CHKERRQ(ierr);
40086d74e61SPeter Brune     *changed = PETSC_FALSE;
40186d74e61SPeter Brune     PetscFunctionReturn(0);
40286d74e61SPeter Brune   }
40386d74e61SPeter Brune 
40486d74e61SPeter Brune   ierr = VecDot(Y,Ylast,&dot);CHKERRQ(ierr);
40586d74e61SPeter Brune   ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr);
40686d74e61SPeter Brune   ierr = VecNorm(Ylast,NORM_2,&ylastnorm);CHKERRQ(ierr);
40786d74e61SPeter Brune   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
40886d74e61SPeter Brune   theta = acos((double)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0));
40986d74e61SPeter Brune   angle_radians = angle * PETSC_PI / 180.;
41086d74e61SPeter Brune   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
41186d74e61SPeter Brune     /* Modify the step Y */
41286d74e61SPeter Brune     PetscReal alpha,ydiffnorm;
41386d74e61SPeter Brune     ierr = VecAXPY(Ylast,-1.0,Y);CHKERRQ(ierr);
41486d74e61SPeter Brune     ierr = VecNorm(Ylast,NORM_2,&ydiffnorm);CHKERRQ(ierr);
41586d74e61SPeter Brune     alpha = ylastnorm / ydiffnorm;
41686d74e61SPeter Brune     ierr = VecCopy(Y,Ylast);CHKERRQ(ierr);
41786d74e61SPeter Brune     ierr = VecScale(Y,alpha);CHKERRQ(ierr);
41886d74e61SPeter Brune     ierr = PetscInfo3(snes,"Angle %G degrees less than threshold %G, corrected step by alpha=%G\n",theta*180./PETSC_PI,angle,alpha);CHKERRQ(ierr);
41986d74e61SPeter Brune   } else {
42086d74e61SPeter Brune     ierr = PetscInfo2(snes,"Angle %G degrees exceeds threshold %G, no correction applied\n",theta*180./PETSC_PI,angle);CHKERRQ(ierr);
42186d74e61SPeter Brune     ierr = VecCopy(Y,Ylast);CHKERRQ(ierr);
42286d74e61SPeter Brune     *changed = PETSC_FALSE;
423bf7f4e0aSPeter Brune   }
424bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
425bf7f4e0aSPeter Brune }
426bf7f4e0aSPeter Brune 
427bf7f4e0aSPeter Brune #undef __FUNCT__
428f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchApply"
429f40b411bSPeter Brune /*@
430cd7522eaSPeter Brune    SNESLineSearchApply - Computes the line-search update.
431f40b411bSPeter Brune 
432f1c6b773SPeter Brune    Collective on SNESLineSearch
433f40b411bSPeter Brune 
434f40b411bSPeter Brune    Input Parameters:
435f40b411bSPeter Brune +  linesearch - The linesearch instance.
436f40b411bSPeter Brune .  X - The current solution.
437f40b411bSPeter Brune .  F - The current function.
438f40b411bSPeter Brune .  fnorm - The current norm.
439f40b411bSPeter Brune .  Y - The search direction.
440f40b411bSPeter Brune 
441f40b411bSPeter Brune    Output Parameters:
442f40b411bSPeter Brune +  X - The new solution.
443f40b411bSPeter Brune .  F - The new function.
444f40b411bSPeter Brune -  fnorm - The new function norm.
445f40b411bSPeter Brune 
446cd7522eaSPeter Brune    Options Database Keys:
447cd7522eaSPeter Brune + -snes_linesearch_type - The Line search method
448cd7522eaSPeter Brune . -snes_linesearch_monitor - Print progress of line searches
449cd7522eaSPeter Brune . -snes_linesearch_damping - The linesearch damping parameter.
450cd7522eaSPeter Brune . -snes_linesearch_norms   - Turn on/off the linesearch norms
451cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess.
452cd7522eaSPeter Brune - -snes_linesearch_max_it - The number of iterations for iterative line searches.
453cd7522eaSPeter Brune 
454cd7522eaSPeter Brune    Notes:
455cd7522eaSPeter Brune    This is typically called from within a SNESSolve() implementation in order to
456cd7522eaSPeter Brune    help with convergence of the nonlinear method.  Various SNES types use line searches
457cd7522eaSPeter Brune    in different ways, but the overarching theme is that a line search is used to determine
458cd7522eaSPeter Brune    an optimal damping parameter of a step at each iteration of the method.  Each
459cd7522eaSPeter Brune    application of the line search may invoke SNESComputeFunction several times, and
460cd7522eaSPeter Brune    therefore may be fairly expensive.
461cd7522eaSPeter Brune 
462f40b411bSPeter Brune    Level: Intermediate
463f40b411bSPeter Brune 
464f1c6b773SPeter Brune    .keywords: SNESLineSearch, Create
465f40b411bSPeter Brune 
466cd7522eaSPeter Brune    .seealso: SNESLineSearchCreate(), SNESLineSearchPreCheck(), SNESLineSearchPostCheck(), SNESSolve(), SNESComputeFunction()
467f40b411bSPeter Brune @*/
468f1c6b773SPeter Brune PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y) {
469bf7f4e0aSPeter Brune   PetscErrorCode ierr;
470bf7f4e0aSPeter Brune   PetscFunctionBegin;
471bf7f4e0aSPeter Brune 
472bf7f4e0aSPeter Brune   /* check the pointers */
473f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
474bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
475bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
476bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
477bf7f4e0aSPeter Brune 
478bf7f4e0aSPeter Brune   linesearch->success = PETSC_TRUE;
479bf7f4e0aSPeter Brune 
480bf7f4e0aSPeter Brune   linesearch->vec_sol = X;
481bf7f4e0aSPeter Brune   linesearch->vec_update = Y;
482bf7f4e0aSPeter Brune   linesearch->vec_func = F;
483bf7f4e0aSPeter Brune 
484f1c6b773SPeter Brune   ierr = SNESLineSearchSetUp(linesearch);CHKERRQ(ierr);
485bf7f4e0aSPeter Brune 
486bf7f4e0aSPeter Brune   if (!linesearch->keeplambda)
487bf7f4e0aSPeter Brune     linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */
488bf7f4e0aSPeter Brune 
489bf7f4e0aSPeter Brune   if (fnorm) {
490bf7f4e0aSPeter Brune     linesearch->fnorm = *fnorm;
491bf7f4e0aSPeter Brune   } else {
492bf7f4e0aSPeter Brune     ierr = VecNorm(F, NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
493bf7f4e0aSPeter Brune   }
494bf7f4e0aSPeter Brune 
495f1c6b773SPeter Brune   ierr = PetscLogEventBegin(SNESLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
496bf7f4e0aSPeter Brune 
497bf7f4e0aSPeter Brune   ierr = (*linesearch->ops->apply)(linesearch);CHKERRQ(ierr);
498bf7f4e0aSPeter Brune 
499f1c6b773SPeter Brune   ierr = PetscLogEventEnd(SNESLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
500bf7f4e0aSPeter Brune 
501bf7f4e0aSPeter Brune   if (fnorm)
502bf7f4e0aSPeter Brune     *fnorm = linesearch->fnorm;
503bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
504bf7f4e0aSPeter Brune }
505bf7f4e0aSPeter Brune 
506bf7f4e0aSPeter Brune #undef __FUNCT__
507f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchDestroy"
508f40b411bSPeter Brune /*@
509f1c6b773SPeter Brune    SNESLineSearchDestroy - Destroys the line search instance.
510f40b411bSPeter Brune 
511f1c6b773SPeter Brune    Collective on SNESLineSearch
512f40b411bSPeter Brune 
513f40b411bSPeter Brune    Input Parameters:
514f40b411bSPeter Brune .  linesearch - The linesearch instance.
515f40b411bSPeter Brune 
516f40b411bSPeter Brune    Level: Intermediate
517f40b411bSPeter Brune 
518f1c6b773SPeter Brune    .keywords: SNESLineSearch, Create
519f40b411bSPeter Brune 
520cd7522eaSPeter Brune    .seealso: SNESLineSearchCreate(), SNESLineSearchReset(), SNESDestroy()
521f40b411bSPeter Brune @*/
522f1c6b773SPeter Brune PetscErrorCode SNESLineSearchDestroy(SNESLineSearch * linesearch) {
523bf7f4e0aSPeter Brune   PetscErrorCode ierr;
524bf7f4e0aSPeter Brune   PetscFunctionBegin;
525bf7f4e0aSPeter Brune   if (!*linesearch) PetscFunctionReturn(0);
526f1c6b773SPeter Brune   PetscValidHeaderSpecific((*linesearch),SNESLINESEARCH_CLASSID,1);
527bf7f4e0aSPeter Brune   if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = 0; PetscFunctionReturn(0);}
528bf7f4e0aSPeter Brune   ierr = PetscObjectDepublish((*linesearch));CHKERRQ(ierr);
529f1c6b773SPeter Brune   ierr = SNESLineSearchReset(*linesearch);
530bf7f4e0aSPeter Brune   if ((*linesearch)->ops->destroy) {
531bf7f4e0aSPeter Brune     (*linesearch)->ops->destroy(*linesearch);
532bf7f4e0aSPeter Brune   }
533bf7f4e0aSPeter Brune   ierr = PetscViewerDestroy(&(*linesearch)->monitor);CHKERRQ(ierr);
534e7058c64SPeter Brune   ierr = PetscHeaderDestroy(linesearch);CHKERRQ(ierr);
535bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
536bf7f4e0aSPeter Brune }
537bf7f4e0aSPeter Brune 
538bf7f4e0aSPeter Brune #undef __FUNCT__
539f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetMonitor"
540f40b411bSPeter Brune /*@
541cd7522eaSPeter Brune    SNESLineSearchSetMonitor - Turns on/off printing useful information and debugging output about the line search.
542bf7f4e0aSPeter Brune 
543bf7f4e0aSPeter Brune    Input Parameters:
544bf7f4e0aSPeter Brune +  snes - nonlinear context obtained from SNESCreate()
545bf7f4e0aSPeter Brune -  flg - PETSC_TRUE to monitor the line search
546bf7f4e0aSPeter Brune 
547bf7f4e0aSPeter Brune    Logically Collective on SNES
548bf7f4e0aSPeter Brune 
549bf7f4e0aSPeter Brune    Options Database:
550cd7522eaSPeter Brune .   -snes_linesearch_monitor - enables the monitor.
551bf7f4e0aSPeter Brune 
552bf7f4e0aSPeter Brune    Level: intermediate
553bf7f4e0aSPeter Brune 
554bf7f4e0aSPeter Brune 
555cd7522eaSPeter Brune .seealso: SNESLineSearchGetMonitor(), PetscViewer
556bf7f4e0aSPeter Brune @*/
557f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetMonitor(SNESLineSearch linesearch, PetscBool flg)
558bf7f4e0aSPeter Brune {
559bf7f4e0aSPeter Brune 
560bf7f4e0aSPeter Brune   PetscErrorCode ierr;
561bf7f4e0aSPeter Brune   PetscFunctionBegin;
562bf7f4e0aSPeter Brune   if (flg && !linesearch->monitor) {
563bf7f4e0aSPeter Brune     ierr = PetscViewerASCIIOpen(((PetscObject)linesearch)->comm,"stdout",&linesearch->monitor);CHKERRQ(ierr);
564bf7f4e0aSPeter Brune   } else if (!flg && linesearch->monitor) {
565bf7f4e0aSPeter Brune     ierr = PetscViewerDestroy(&linesearch->monitor);CHKERRQ(ierr);
566bf7f4e0aSPeter Brune   }
567bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
568bf7f4e0aSPeter Brune }
569bf7f4e0aSPeter Brune 
570bf7f4e0aSPeter Brune #undef __FUNCT__
571f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetMonitor"
572f40b411bSPeter Brune /*@
573cd7522eaSPeter Brune    SNESLineSearchGetMonitor - Gets the PetscViewer instance for the line search monitor.
5746a388c36SPeter Brune 
575f40b411bSPeter Brune    Input Parameters:
576f40b411bSPeter Brune .  linesearch - linesearch context.
577f40b411bSPeter Brune 
578f40b411bSPeter Brune    Input Parameters:
579f40b411bSPeter Brune .  monitor - monitor context.
580f40b411bSPeter Brune 
581f40b411bSPeter Brune    Logically Collective on SNES
582f40b411bSPeter Brune 
583f40b411bSPeter Brune 
584f40b411bSPeter Brune    Options Database Keys:
585cd7522eaSPeter Brune .   -snes_linesearch_monitor - enables the monitor.
586f40b411bSPeter Brune 
587f40b411bSPeter Brune    Level: intermediate
588f40b411bSPeter Brune 
589f40b411bSPeter Brune 
590cd7522eaSPeter Brune .seealso: SNESLineSearchSetMonitor(), PetscViewer
591f40b411bSPeter Brune @*/
592f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetMonitor(SNESLineSearch linesearch, PetscViewer *monitor)
5936a388c36SPeter Brune {
5946a388c36SPeter Brune 
5956a388c36SPeter Brune   PetscFunctionBegin;
596f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
5976a388c36SPeter Brune   if (monitor) {
5986a388c36SPeter Brune     PetscValidPointer(monitor, 2);
5996a388c36SPeter Brune     *monitor = linesearch->monitor;
6006a388c36SPeter Brune   }
6016a388c36SPeter Brune   PetscFunctionReturn(0);
6026a388c36SPeter Brune }
6036a388c36SPeter Brune 
6046a388c36SPeter Brune #undef __FUNCT__
605f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetFromOptions"
606f40b411bSPeter Brune /*@
607f1c6b773SPeter Brune    SNESLineSearchSetFromOptions - Sets options for the line search
608f40b411bSPeter Brune 
609f40b411bSPeter Brune    Input Parameters:
610f40b411bSPeter Brune .  linesearch - linesearch context.
611f40b411bSPeter Brune 
612f40b411bSPeter Brune    Options Database Keys:
613cd7522eaSPeter Brune + -snes_linesearch_type - The Line search method
614cd7522eaSPeter Brune . -snes_linesearch_monitor - Print progress of line searches
615cd7522eaSPeter Brune . -snes_linesearch_damping - The linesearch damping parameter.
616cd7522eaSPeter Brune . -snes_linesearch_norms   - Turn on/off the linesearch norms
617cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess.
61859405d5eSPeter Brune . -snes_linesearch_order - The polynomial order of approximation used in the linesearch
619cd7522eaSPeter Brune - -snes_linesearch_max_it - The number of iterations for iterative line searches.
620f40b411bSPeter Brune 
62159405d5eSPeter Brune 
622f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
623f40b411bSPeter Brune 
624f40b411bSPeter Brune    Level: intermediate
625f40b411bSPeter Brune 
626f40b411bSPeter Brune 
627f1c6b773SPeter Brune .seealso: SNESLineSearchCreate()
628f40b411bSPeter Brune @*/
629f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch) {
630bf7f4e0aSPeter Brune   PetscErrorCode ierr;
63159405d5eSPeter Brune   const char       *orders[] = {"linear", "quadratic", "cubic"};
63259405d5eSPeter Brune   PetscInt         indx = 0;
633f1c6b773SPeter Brune   const char     *deft = SNES_LINESEARCH_BASIC;
634bf7f4e0aSPeter Brune   char           type[256];
635bf7f4e0aSPeter Brune   PetscBool      flg, set;
636bf7f4e0aSPeter Brune   PetscFunctionBegin;
637f1c6b773SPeter Brune   if (!SNESLineSearchRegisterAllCalled) {ierr = SNESLineSearchRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
638bf7f4e0aSPeter Brune 
639bf7f4e0aSPeter Brune   ierr = PetscObjectOptionsBegin((PetscObject)linesearch);CHKERRQ(ierr);
640bf7f4e0aSPeter Brune   if (((PetscObject)linesearch)->type_name) {
641bf7f4e0aSPeter Brune     deft = ((PetscObject)linesearch)->type_name;
642bf7f4e0aSPeter Brune   }
6437a35526eSPeter Brune   ierr = PetscOptionsList("-snes_linesearch_type","Line-search method","SNESLineSearchSetType",SNESLineSearchList,deft,type,256,&flg);CHKERRQ(ierr);
644bf7f4e0aSPeter Brune   if (flg) {
645f1c6b773SPeter Brune     ierr = SNESLineSearchSetType(linesearch,type);CHKERRQ(ierr);
646bf7f4e0aSPeter Brune   } else if (!((PetscObject)linesearch)->type_name) {
647f1c6b773SPeter Brune     ierr = SNESLineSearchSetType(linesearch,deft);CHKERRQ(ierr);
648bf7f4e0aSPeter Brune   }
649bf7f4e0aSPeter Brune   if (linesearch->ops->setfromoptions) {
650bf7f4e0aSPeter Brune     (*linesearch->ops->setfromoptions)(linesearch);CHKERRQ(ierr);
651bf7f4e0aSPeter Brune   }
652bf7f4e0aSPeter Brune 
6537a35526eSPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_monitor","Print progress of line searches","SNESSNESLineSearchSetMonitor",
654bf7f4e0aSPeter Brune                           linesearch->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
655f1c6b773SPeter Brune   if (set) {ierr = SNESLineSearchSetMonitor(linesearch,flg);CHKERRQ(ierr);}
656bf7f4e0aSPeter Brune 
6577a35526eSPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_damping","Line search damping and initial step guess","SNESLineSearchSetDamping",linesearch->damping,&linesearch->damping,0);CHKERRQ(ierr);
6587a35526eSPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_rtol","Tolerance for iterative line search","SNESLineSearchSetRTolerance",linesearch->rtol,&linesearch->rtol,0);CHKERRQ(ierr);
6597a35526eSPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_norms","Compute final norms in line search","SNESLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,0);CHKERRQ(ierr);
6607a35526eSPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_keeplambda","Use previous lambda as damping","SNESLineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,0);CHKERRQ(ierr);
6617a35526eSPeter Brune   ierr = PetscOptionsInt("-snes_linesearch_max_it","Maximum iterations for iterative line searches","",linesearch->max_its,&linesearch->max_its,0);CHKERRQ(ierr);
6627a35526eSPeter Brune 
6637a35526eSPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_precheck_picard","Use a correction that sometimes improves convergence of Picard iteration","SNESLineSearchPreCheckPicard",flg,&flg,&set);CHKERRQ(ierr);
6647a35526eSPeter Brune   if (set) {
6657a35526eSPeter Brune     if (flg) {
6667a35526eSPeter Brune       linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */
6677a35526eSPeter Brune       ierr = PetscOptionsReal("-snes_linesearch_precheck_picard_angle","Maximum angle at which to activate the correction",
6687a35526eSPeter Brune                               "none",linesearch->precheck_picard_angle,&linesearch->precheck_picard_angle,PETSC_NULL);CHKERRQ(ierr);
6697a35526eSPeter Brune       ierr = SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&linesearch->precheck_picard_angle);CHKERRQ(ierr);
6707a35526eSPeter Brune     } else {
6717a35526eSPeter Brune       ierr = SNESLineSearchSetPreCheck(linesearch,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
6727a35526eSPeter Brune     }
6737a35526eSPeter Brune   }
67459405d5eSPeter Brune   ierr = PetscOptionsEList("-snes_linesearch_order",  "Order of approximation", "SNESLineSearch", orders,3,orders[(PetscInt)linesearch->order],&indx,&flg);CHKERRQ(ierr);
67559405d5eSPeter Brune   if (flg) {
67659405d5eSPeter Brune     switch (indx) {
67759405d5eSPeter Brune     case 0: linesearch->order = SNES_LINESEARCH_LINEAR;
67859405d5eSPeter Brune       break;
67959405d5eSPeter Brune     case 1: linesearch->order = SNES_LINESEARCH_QUADRATIC;
68059405d5eSPeter Brune       break;
68159405d5eSPeter Brune     case 2: linesearch->order = SNES_LINESEARCH_CUBIC;
68259405d5eSPeter Brune       break;
68359405d5eSPeter Brune     }
68459405d5eSPeter Brune   }
68559405d5eSPeter Brune 
6867a35526eSPeter Brune 
687bf7f4e0aSPeter Brune   ierr = PetscObjectProcessOptionsHandlers((PetscObject)linesearch);CHKERRQ(ierr);
688bf7f4e0aSPeter Brune   ierr = PetscOptionsEnd();CHKERRQ(ierr);
689bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
690bf7f4e0aSPeter Brune }
691bf7f4e0aSPeter Brune 
692bf7f4e0aSPeter Brune #undef __FUNCT__
693f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchView"
694f40b411bSPeter Brune /*@
695cd7522eaSPeter Brune    SNESLineSearchView - Prints useful information about the line search not
696cd7522eaSPeter Brune    related to an individual call.
697f40b411bSPeter Brune 
698f40b411bSPeter Brune    Input Parameters:
699f40b411bSPeter Brune .  linesearch - linesearch context.
700f40b411bSPeter Brune 
701f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
702f40b411bSPeter Brune 
703f40b411bSPeter Brune    Level: intermediate
704f40b411bSPeter Brune 
705f1c6b773SPeter Brune .seealso: SNESLineSearchCreate()
706f40b411bSPeter Brune @*/
7077f1410a3SPeter Brune PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer) {
7087f1410a3SPeter Brune   PetscErrorCode ierr;
7097f1410a3SPeter Brune   PetscBool      iascii;
710bf7f4e0aSPeter Brune   PetscFunctionBegin;
7117f1410a3SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
7127f1410a3SPeter Brune   if (!viewer) {
7137f1410a3SPeter Brune     ierr = PetscViewerASCIIGetStdout(((PetscObject)linesearch)->comm,&viewer);CHKERRQ(ierr);
7147f1410a3SPeter Brune   }
7157f1410a3SPeter Brune   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
7167f1410a3SPeter Brune   PetscCheckSameComm(linesearch,1,viewer,2);
717f40b411bSPeter Brune 
7187f1410a3SPeter Brune   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
7197f1410a3SPeter Brune   if (iascii) {
7207f1410a3SPeter Brune     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)linesearch,viewer,"SNESLineSearch Object");CHKERRQ(ierr);
7217f1410a3SPeter Brune     if (linesearch->ops->view) {
7227f1410a3SPeter Brune       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
7237f1410a3SPeter Brune       ierr = (*linesearch->ops->view)(linesearch,viewer);CHKERRQ(ierr);
7247f1410a3SPeter Brune       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
7257f1410a3SPeter Brune     }
7267f1410a3SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  maxstep=%G, minlambda=%G\n", linesearch->maxstep,linesearch->steptol);CHKERRQ(ierr);
7277f1410a3SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%G, absolute=%G, lambda=%G\n", linesearch->rtol,linesearch->atol,linesearch->ltol);CHKERRQ(ierr);
7287f1410a3SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D\n", linesearch->max_its);CHKERRQ(ierr);
7297f1410a3SPeter Brune     if (linesearch->ops->precheckstep) {
7307f1410a3SPeter Brune       if (linesearch->ops->precheckstep == SNESLineSearchPreCheckPicard) {
7317f1410a3SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  using precheck step to speed up Picard convergence\n", linesearch->max_its);CHKERRQ(ierr);
7327f1410a3SPeter Brune       } else {
7337f1410a3SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  using user-defined precheck step\n", linesearch->max_its);CHKERRQ(ierr);
7347f1410a3SPeter Brune       }
7357f1410a3SPeter Brune     }
7367f1410a3SPeter Brune     if (linesearch->ops->postcheckstep) {
7377f1410a3SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer,"  using user-defined postcheck step\n", linesearch->max_its);CHKERRQ(ierr);
7387f1410a3SPeter Brune     }
7397f1410a3SPeter Brune   }
740bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
741bf7f4e0aSPeter Brune }
742bf7f4e0aSPeter Brune 
743bf7f4e0aSPeter Brune #undef __FUNCT__
744f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetType"
745ea5d4fccSPeter Brune /*@C
746f1c6b773SPeter Brune    SNESLineSearchSetType - Sets the linesearch type
747f40b411bSPeter Brune 
748f40b411bSPeter Brune    Input Parameters:
749f40b411bSPeter Brune +  linesearch - linesearch context.
750f40b411bSPeter Brune -  type - The type of line search to be used
751f40b411bSPeter Brune 
752cd7522eaSPeter Brune    Available Types:
753cd7522eaSPeter Brune +  basic - Simple damping line search.
754cd7522eaSPeter Brune .  bt - Backtracking line search over the L2 norm of the function.
755cd7522eaSPeter Brune .  l2 - Secant line search over the L2 norm of the function.
756cd7522eaSPeter Brune .  cp - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x).
757cd7522eaSPeter Brune -  shell - User provided SNESLineSearch implementation.
758cd7522eaSPeter Brune 
759f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
760f40b411bSPeter Brune 
761f40b411bSPeter Brune    Level: intermediate
762f40b411bSPeter Brune 
763f40b411bSPeter Brune 
764f1c6b773SPeter Brune .seealso: SNESLineSearchCreate()
765f40b411bSPeter Brune @*/
766f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, const SNESLineSearchType type)
767bf7f4e0aSPeter Brune {
768bf7f4e0aSPeter Brune 
769f1c6b773SPeter Brune   PetscErrorCode ierr,(*r)(SNESLineSearch);
770bf7f4e0aSPeter Brune   PetscBool      match;
771bf7f4e0aSPeter Brune 
772bf7f4e0aSPeter Brune   PetscFunctionBegin;
773f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
774bf7f4e0aSPeter Brune   PetscValidCharPointer(type,2);
775bf7f4e0aSPeter Brune 
776bf7f4e0aSPeter Brune   ierr = PetscTypeCompare((PetscObject)linesearch,type,&match);CHKERRQ(ierr);
777bf7f4e0aSPeter Brune   if (match) PetscFunctionReturn(0);
778bf7f4e0aSPeter Brune 
779f1c6b773SPeter Brune   ierr =  PetscFListFind(SNESLineSearchList,((PetscObject)linesearch)->comm,type,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr);
780bf7f4e0aSPeter Brune   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type);
781bf7f4e0aSPeter Brune   /* Destroy the previous private linesearch context */
782bf7f4e0aSPeter Brune   if (linesearch->ops->destroy) {
783bf7f4e0aSPeter Brune     ierr = (*(linesearch)->ops->destroy)(linesearch);CHKERRQ(ierr);
784bf7f4e0aSPeter Brune     linesearch->ops->destroy = PETSC_NULL;
785bf7f4e0aSPeter Brune   }
786f1c6b773SPeter Brune   /* Reinitialize function pointers in SNESLineSearchOps structure */
787bf7f4e0aSPeter Brune   linesearch->ops->apply          = 0;
788bf7f4e0aSPeter Brune   linesearch->ops->view           = 0;
789bf7f4e0aSPeter Brune   linesearch->ops->setfromoptions = 0;
790bf7f4e0aSPeter Brune   linesearch->ops->destroy        = 0;
791bf7f4e0aSPeter Brune 
792bf7f4e0aSPeter Brune   ierr = PetscObjectChangeTypeName((PetscObject)linesearch,type);CHKERRQ(ierr);
793bf7f4e0aSPeter Brune   ierr = (*r)(linesearch);CHKERRQ(ierr);
794bf7f4e0aSPeter Brune #if defined(PETSC_HAVE_AMS)
795bf7f4e0aSPeter Brune   if (PetscAMSPublishAll) {
796bf7f4e0aSPeter Brune     ierr = PetscObjectAMSPublish((PetscObject)linesearch);CHKERRQ(ierr);
797bf7f4e0aSPeter Brune   }
798bf7f4e0aSPeter Brune #endif
799bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
800bf7f4e0aSPeter Brune }
801bf7f4e0aSPeter Brune 
802bf7f4e0aSPeter Brune #undef __FUNCT__
803f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetSNES"
804f40b411bSPeter Brune /*@
805f1c6b773SPeter Brune    SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation
806f40b411bSPeter Brune 
807f40b411bSPeter Brune    Input Parameters:
808f40b411bSPeter Brune +  linesearch - linesearch context.
809f40b411bSPeter Brune -  snes - The snes instance
810f40b411bSPeter Brune 
811f40b411bSPeter Brune    Level: intermediate
812f40b411bSPeter Brune 
813f40b411bSPeter Brune 
814cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
815f40b411bSPeter Brune @*/
816f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes){
817bf7f4e0aSPeter Brune   PetscFunctionBegin;
818f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
819bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,2);
820bf7f4e0aSPeter Brune   linesearch->snes = snes;
821bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
822bf7f4e0aSPeter Brune }
823bf7f4e0aSPeter Brune 
824bf7f4e0aSPeter Brune #undef __FUNCT__
825f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetSNES"
826f40b411bSPeter Brune /*@
827f1c6b773SPeter Brune    SNESLineSearchGetSNES - Gets the SNES for the linesearch for function evaluation
828f40b411bSPeter Brune 
829f40b411bSPeter Brune    Input Parameters:
830f40b411bSPeter Brune .  linesearch - linesearch context.
831f40b411bSPeter Brune 
832f40b411bSPeter Brune    Output Parameters:
833f40b411bSPeter Brune .  snes - The snes instance
834f40b411bSPeter Brune 
835f40b411bSPeter Brune    Level: intermediate
836f40b411bSPeter Brune 
837cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
838f40b411bSPeter Brune @*/
839f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes){
840bf7f4e0aSPeter Brune   PetscFunctionBegin;
841f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
8426a388c36SPeter Brune   PetscValidPointer(snes, 2);
843bf7f4e0aSPeter Brune   *snes = linesearch->snes;
844bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
845bf7f4e0aSPeter Brune }
846bf7f4e0aSPeter Brune 
8476a388c36SPeter Brune #undef __FUNCT__
848f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetLambda"
849f40b411bSPeter Brune /*@
850f1c6b773SPeter Brune    SNESLineSearchGetLambda - Gets the last linesearch steplength discovered.
851f40b411bSPeter Brune 
852f40b411bSPeter Brune    Input Parameters:
853f40b411bSPeter Brune .  linesearch - linesearch context.
854f40b411bSPeter Brune 
855f40b411bSPeter Brune    Output Parameters:
856cd7522eaSPeter Brune .  lambda - The last steplength computed during SNESLineSearchApply()
857f40b411bSPeter Brune 
858f40b411bSPeter Brune    Level: intermediate
859f40b411bSPeter Brune 
860cd7522eaSPeter Brune .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetDamping(), SNESLineSearchApply()
861f40b411bSPeter Brune @*/
862f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda)
8636a388c36SPeter Brune {
8646a388c36SPeter Brune   PetscFunctionBegin;
865f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
8666a388c36SPeter Brune   PetscValidPointer(lambda, 2);
8676a388c36SPeter Brune   *lambda = linesearch->lambda;
8686a388c36SPeter Brune   PetscFunctionReturn(0);
8696a388c36SPeter Brune }
8706a388c36SPeter Brune 
8716a388c36SPeter Brune #undef __FUNCT__
872f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetLambda"
873f40b411bSPeter Brune /*@
874f1c6b773SPeter Brune    SNESLineSearchSetLambda - Sets the linesearch steplength.
875f40b411bSPeter Brune 
876f40b411bSPeter Brune    Input Parameters:
877f40b411bSPeter Brune +  linesearch - linesearch context.
878f40b411bSPeter Brune -  lambda - The last steplength.
879f40b411bSPeter Brune 
880cd7522eaSPeter Brune    Notes:
881cd7522eaSPeter Brune    This routine is typically used within implementations of SNESLineSearchApply
882cd7522eaSPeter Brune    to set the final steplength.  This routine (and SNESLineSearchGetLambda()) were
883cd7522eaSPeter Brune    added in order to facilitate Quasi-Newton methods that use the previous steplength
884cd7522eaSPeter Brune    as an inner scaling parameter.
885cd7522eaSPeter Brune 
886f40b411bSPeter Brune    Level: intermediate
887f40b411bSPeter Brune 
888f1c6b773SPeter Brune .seealso: SNESLineSearchGetLambda()
889f40b411bSPeter Brune @*/
890f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
8916a388c36SPeter Brune {
8926a388c36SPeter Brune   PetscFunctionBegin;
893f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
8946a388c36SPeter Brune   linesearch->lambda = lambda;
8956a388c36SPeter Brune   PetscFunctionReturn(0);
8966a388c36SPeter Brune }
8976a388c36SPeter Brune 
8986a388c36SPeter Brune #undef  __FUNCT__
899f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetTolerances"
900f40b411bSPeter Brune /*@
901f1c6b773SPeter Brune    SNESLineSearchGetTolerances - Gets the tolerances for the method
902f40b411bSPeter Brune 
903f40b411bSPeter Brune    Input Parameters:
904f40b411bSPeter Brune .  linesearch - linesearch context.
905f40b411bSPeter Brune 
906f40b411bSPeter Brune    Output Parameters:
907516fe3c3SPeter Brune +  steptol - The minimum steplength
9086cc8e53bSPeter Brune .  maxstep - The maximum steplength
909516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
910516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
911516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
912516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
913f40b411bSPeter Brune 
914516fe3c3SPeter Brune    Level: advanced
915516fe3c3SPeter Brune 
916f1c6b773SPeter Brune .seealso: SNESLineSearchSetTolerances()
917f40b411bSPeter Brune @*/
918f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
9196a388c36SPeter Brune {
9206a388c36SPeter Brune   PetscFunctionBegin;
921f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
922516fe3c3SPeter Brune   if (steptol) {
9236a388c36SPeter Brune     PetscValidPointer(steptol, 2);
9246a388c36SPeter Brune     *steptol = linesearch->steptol;
925516fe3c3SPeter Brune   }
926516fe3c3SPeter Brune   if (maxstep) {
927516fe3c3SPeter Brune     PetscValidPointer(maxstep, 3);
928516fe3c3SPeter Brune     *maxstep = linesearch->maxstep;
929516fe3c3SPeter Brune   }
930516fe3c3SPeter Brune   if (rtol) {
931516fe3c3SPeter Brune     PetscValidPointer(rtol, 4);
932516fe3c3SPeter Brune     *rtol = linesearch->rtol;
933516fe3c3SPeter Brune   }
934516fe3c3SPeter Brune   if (atol) {
935516fe3c3SPeter Brune     PetscValidPointer(atol, 5);
936516fe3c3SPeter Brune     *atol = linesearch->atol;
937516fe3c3SPeter Brune   }
938516fe3c3SPeter Brune   if (ltol) {
939516fe3c3SPeter Brune     PetscValidPointer(ltol, 6);
940516fe3c3SPeter Brune     *ltol = linesearch->ltol;
941516fe3c3SPeter Brune   }
942516fe3c3SPeter Brune   if (max_its) {
943516fe3c3SPeter Brune     PetscValidPointer(max_its, 7);
944516fe3c3SPeter Brune     *max_its = linesearch->max_its;
945516fe3c3SPeter Brune   }
9466a388c36SPeter Brune   PetscFunctionReturn(0);
9476a388c36SPeter Brune }
9486a388c36SPeter Brune 
9496a388c36SPeter Brune #undef  __FUNCT__
950f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetTolerances"
951f40b411bSPeter Brune /*@
952f1c6b773SPeter Brune    SNESLineSearchSetTolerances - Sets the tolerances for the method
953f40b411bSPeter Brune 
954f40b411bSPeter Brune    Input Parameters:
955516fe3c3SPeter Brune +  linesearch - linesearch context.
956516fe3c3SPeter Brune .  steptol - The minimum steplength
9576cc8e53bSPeter Brune .  maxstep - The maximum steplength
958516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
959516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
960516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
961516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
962f40b411bSPeter Brune 
963f40b411bSPeter Brune 
964516fe3c3SPeter Brune    Level: advanced
965516fe3c3SPeter Brune 
966f1c6b773SPeter Brune .seealso: SNESLineSearchGetTolerances()
967f40b411bSPeter Brune @*/
968f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its)
9696a388c36SPeter Brune {
9706a388c36SPeter Brune   PetscFunctionBegin;
971f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
972d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,steptol,2);
973d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,maxstep,3);
974d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,rtol,4);
975d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,atol,5);
976d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,ltol,6);
977d3952184SSatish Balay   PetscValidLogicalCollectiveInt(linesearch,max_its,7);
978d3952184SSatish Balay 
979d3952184SSatish Balay   if ( steptol!= PETSC_DEFAULT) {
980d3952184SSatish Balay     if ( steptol < 0.0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Minimum step length %G must be non-negative",steptol);
9816a388c36SPeter Brune     linesearch->steptol = steptol;
982d3952184SSatish Balay   }
983d3952184SSatish Balay 
984d3952184SSatish Balay   if ( maxstep!= PETSC_DEFAULT) {
985d3952184SSatish Balay     if ( maxstep < 0.0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Maximum step length %G must be non-negative",maxstep);
986516fe3c3SPeter Brune     linesearch->maxstep = maxstep;
987d3952184SSatish Balay   }
988d3952184SSatish Balay 
989d3952184SSatish Balay   if (rtol != PETSC_DEFAULT) {
990d3952184SSatish 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);
991516fe3c3SPeter Brune     linesearch->rtol = rtol;
992d3952184SSatish Balay   }
993d3952184SSatish Balay 
994d3952184SSatish Balay   if (atol != PETSC_DEFAULT) {
995d3952184SSatish Balay     if (atol < 0.0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %G must be non-negative",atol);
996516fe3c3SPeter Brune     linesearch->atol = atol;
997d3952184SSatish Balay   }
998d3952184SSatish Balay 
999d3952184SSatish Balay   if (ltol != PETSC_DEFAULT) {
1000d3952184SSatish Balay     if (ltol < 0.0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Labmda tolerance %G must be non-negative",ltol);
1001516fe3c3SPeter Brune   linesearch->ltol = ltol;
1002d3952184SSatish Balay   }
1003d3952184SSatish Balay 
1004d3952184SSatish Balay   if (max_its != PETSC_DEFAULT) {
1005d3952184SSatish Balay     if (max_its < 0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",max_its);
1006516fe3c3SPeter Brune     linesearch->max_its = max_its;
1007d3952184SSatish Balay   }
1008d3952184SSatish Balay 
10096a388c36SPeter Brune   PetscFunctionReturn(0);
10106a388c36SPeter Brune }
10116a388c36SPeter Brune 
1012516fe3c3SPeter Brune 
10136a388c36SPeter Brune #undef __FUNCT__
1014f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetDamping"
1015f40b411bSPeter Brune /*@
1016f1c6b773SPeter Brune    SNESLineSearchGetDamping - Gets the line search damping parameter.
1017f40b411bSPeter Brune 
1018f40b411bSPeter Brune    Input Parameters:
1019f40b411bSPeter Brune .  linesearch - linesearch context.
1020f40b411bSPeter Brune 
1021f40b411bSPeter Brune    Output Parameters:
1022f40b411bSPeter Brune .  damping - The damping parameter.
1023f40b411bSPeter Brune 
1024f40b411bSPeter Brune    Level: intermediate
1025f40b411bSPeter Brune 
1026f1c6b773SPeter Brune .seealso: SNESLineSearchGetStepTolerance()
1027f40b411bSPeter Brune @*/
1028f40b411bSPeter Brune 
1029f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping)
10306a388c36SPeter Brune {
10316a388c36SPeter Brune   PetscFunctionBegin;
1032f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
10336a388c36SPeter Brune   PetscValidPointer(damping, 2);
10346a388c36SPeter Brune   *damping = linesearch->damping;
10356a388c36SPeter Brune   PetscFunctionReturn(0);
10366a388c36SPeter Brune }
10376a388c36SPeter Brune 
10386a388c36SPeter Brune #undef __FUNCT__
1039f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetDamping"
1040f40b411bSPeter Brune /*@
1041f1c6b773SPeter Brune    SNESLineSearchSetDamping - Sets the line search damping paramter.
1042f40b411bSPeter Brune 
1043f40b411bSPeter Brune    Input Parameters:
1044f40b411bSPeter Brune .  linesearch - linesearch context.
1045f40b411bSPeter Brune .  damping - The damping parameter.
1046f40b411bSPeter Brune 
1047f40b411bSPeter Brune    Level: intermediate
1048f40b411bSPeter Brune 
1049cd7522eaSPeter Brune    Notes:
1050cd7522eaSPeter Brune    The basic line search merely takes the update step scaled by the damping parameter.
1051cd7522eaSPeter Brune    The use of the damping parameter in the l2 and cp line searches is much more subtle;
1052cd7522eaSPeter Brune    it is used as a starting point in calculating the secant step; however, the eventual
1053cd7522eaSPeter Brune    step may be of greater length than the damping parameter.  In the bt line search it is
1054cd7522eaSPeter Brune    used as the maximum possible step length, as the bt line search only backtracks.
1055cd7522eaSPeter Brune 
1056f1c6b773SPeter Brune .seealso: SNESLineSearchGetDamping()
1057f40b411bSPeter Brune @*/
1058f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping)
10596a388c36SPeter Brune {
10606a388c36SPeter Brune   PetscFunctionBegin;
1061f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
10626a388c36SPeter Brune   linesearch->damping = damping;
10636a388c36SPeter Brune   PetscFunctionReturn(0);
10646a388c36SPeter Brune }
10656a388c36SPeter Brune 
10666a388c36SPeter Brune #undef __FUNCT__
106759405d5eSPeter Brune #define __FUNCT__ "SNESLineSearchGetOrder"
106859405d5eSPeter Brune /*@
106959405d5eSPeter Brune    SNESLineSearchGetOrder - Gets the line search approximation order.
107059405d5eSPeter Brune 
107159405d5eSPeter Brune    Input Parameters:
107259405d5eSPeter Brune .  linesearch - linesearch context.
107359405d5eSPeter Brune 
107459405d5eSPeter Brune    Output Parameters:
107559405d5eSPeter Brune .  order - The order.
107659405d5eSPeter Brune 
107759405d5eSPeter Brune    Level: intermediate
107859405d5eSPeter Brune 
107959405d5eSPeter Brune .seealso: SNESLineSearchSetOrder()
108059405d5eSPeter Brune @*/
108159405d5eSPeter Brune 
108259405d5eSPeter Brune PetscErrorCode  SNESLineSearchGetOrder(SNESLineSearch linesearch,SNESLineSearchOrder *order)
108359405d5eSPeter Brune {
108459405d5eSPeter Brune   PetscFunctionBegin;
108559405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
108659405d5eSPeter Brune   PetscValidPointer(order, 2);
108759405d5eSPeter Brune   *order = linesearch->order;
108859405d5eSPeter Brune   PetscFunctionReturn(0);
108959405d5eSPeter Brune }
109059405d5eSPeter Brune 
109159405d5eSPeter Brune #undef __FUNCT__
109259405d5eSPeter Brune #define __FUNCT__ "SNESLineSearchSetOrder"
109359405d5eSPeter Brune /*@
109459405d5eSPeter Brune    SNESLineSearchSetOrder - Sets the line search damping paramter.
109559405d5eSPeter Brune 
109659405d5eSPeter Brune    Input Parameters:
109759405d5eSPeter Brune .  linesearch - linesearch context.
109859405d5eSPeter Brune .  order - The damping parameter.
109959405d5eSPeter Brune 
110059405d5eSPeter Brune    Level: intermediate
110159405d5eSPeter Brune 
110259405d5eSPeter Brune    Notes:
110359405d5eSPeter Brune    Variable orders are supported by the following line searches:
110459405d5eSPeter Brune .  bt - cubic and quadratic
110559405d5eSPeter Brune 
110659405d5eSPeter Brune .seealso: SNESLineSearchGetOrder()
110759405d5eSPeter Brune @*/
110859405d5eSPeter Brune PetscErrorCode  SNESLineSearchSetOrder(SNESLineSearch linesearch,SNESLineSearchOrder order)
110959405d5eSPeter Brune {
111059405d5eSPeter Brune   PetscFunctionBegin;
111159405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
111259405d5eSPeter Brune   linesearch->order = order;
111359405d5eSPeter Brune   PetscFunctionReturn(0);
111459405d5eSPeter Brune }
111559405d5eSPeter Brune 
111659405d5eSPeter Brune #undef __FUNCT__
1117f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetNorms"
1118f40b411bSPeter Brune /*@
1119f1c6b773SPeter Brune    SNESLineSearchGetNorms - Gets the norms for for X, Y, and F.
1120f40b411bSPeter Brune 
1121f40b411bSPeter Brune    Input Parameters:
1122f40b411bSPeter Brune .  linesearch - linesearch context.
1123f40b411bSPeter Brune 
1124f40b411bSPeter Brune    Output Parameters:
1125f40b411bSPeter Brune +  xnorm - The norm of the current solution
1126f40b411bSPeter Brune .  fnorm - The norm of the current function
1127f40b411bSPeter Brune -  ynorm - The norm of the current update
1128f40b411bSPeter Brune 
1129cd7522eaSPeter Brune    Notes:
1130cd7522eaSPeter Brune    This function is mainly called from SNES implementations.
1131cd7522eaSPeter Brune 
1132f40b411bSPeter Brune    Level: intermediate
1133f40b411bSPeter Brune 
1134f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms() SNESLineSearchGetVecs()
1135f40b411bSPeter Brune @*/
1136f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm)
1137bf7f4e0aSPeter Brune {
1138bf7f4e0aSPeter Brune   PetscFunctionBegin;
1139f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1140bf7f4e0aSPeter Brune   if (xnorm) {
1141bf7f4e0aSPeter Brune     *xnorm = linesearch->xnorm;
1142bf7f4e0aSPeter Brune   }
1143bf7f4e0aSPeter Brune   if (fnorm) {
1144bf7f4e0aSPeter Brune     *fnorm = linesearch->fnorm;
1145bf7f4e0aSPeter Brune   }
1146bf7f4e0aSPeter Brune   if (ynorm) {
1147bf7f4e0aSPeter Brune     *ynorm = linesearch->ynorm;
1148bf7f4e0aSPeter Brune   }
1149bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1150bf7f4e0aSPeter Brune }
1151bf7f4e0aSPeter Brune 
1152e7058c64SPeter Brune #undef __FUNCT__
1153f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetNorms"
1154f40b411bSPeter Brune /*@
1155f1c6b773SPeter Brune    SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F.
1156f40b411bSPeter Brune 
1157f40b411bSPeter Brune    Input Parameters:
1158f40b411bSPeter Brune +  linesearch - linesearch context.
1159f40b411bSPeter Brune .  xnorm - The norm of the current solution
1160f40b411bSPeter Brune .  fnorm - The norm of the current function
1161f40b411bSPeter Brune -  ynorm - The norm of the current update
1162f40b411bSPeter Brune 
1163f40b411bSPeter Brune    Level: intermediate
1164f40b411bSPeter Brune 
1165f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1166f40b411bSPeter Brune @*/
1167f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
11686a388c36SPeter Brune {
11696a388c36SPeter Brune   PetscFunctionBegin;
1170f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
11716a388c36SPeter Brune   linesearch->xnorm = xnorm;
11726a388c36SPeter Brune   linesearch->fnorm = fnorm;
11736a388c36SPeter Brune   linesearch->ynorm = ynorm;
11746a388c36SPeter Brune   PetscFunctionReturn(0);
11756a388c36SPeter Brune }
11766a388c36SPeter Brune 
11776a388c36SPeter Brune #undef __FUNCT__
1178f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchComputeNorms"
1179f40b411bSPeter Brune /*@
1180f1c6b773SPeter Brune    SNESLineSearchComputeNorms - Computes the norms of X, F, and Y.
1181f40b411bSPeter Brune 
1182f40b411bSPeter Brune    Input Parameters:
1183f40b411bSPeter Brune .  linesearch - linesearch context.
1184f40b411bSPeter Brune 
1185f40b411bSPeter Brune    Options Database Keys:
1186cd7522eaSPeter Brune .   -snes_linesearch_norms - turn norm computation on or off.
1187f40b411bSPeter Brune 
1188f40b411bSPeter Brune    Level: intermediate
1189f40b411bSPeter Brune 
1190f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms, SNESLineSearchSetNorms()
1191f40b411bSPeter Brune @*/
1192f1c6b773SPeter Brune PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
11936a388c36SPeter Brune {
11946a388c36SPeter Brune   PetscErrorCode ierr;
11959bd66eb0SPeter Brune   SNES snes;
11966a388c36SPeter Brune   PetscFunctionBegin;
11976a388c36SPeter Brune   if (linesearch->norms) {
11989bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
1199f1c6b773SPeter Brune       ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
12009bd66eb0SPeter Brune       ierr = VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
12019bd66eb0SPeter Brune       ierr = VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
12029bd66eb0SPeter Brune       ierr = (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);CHKERRQ(ierr);
12039bd66eb0SPeter Brune     } else {
12046a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_func,   NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
12056a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_sol,    NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
12066a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
12076a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_func,     NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
12086a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_sol,      NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
12096a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_update,   NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
12106a388c36SPeter Brune     }
12119bd66eb0SPeter Brune   }
12126a388c36SPeter Brune   PetscFunctionReturn(0);
12136a388c36SPeter Brune }
12146a388c36SPeter Brune 
1215*6f263ca3SPeter Brune 
1216*6f263ca3SPeter Brune #undef __FUNCT__
1217*6f263ca3SPeter Brune #define __FUNCT__ "SNESLineSearchSetComputeNorms"
1218*6f263ca3SPeter Brune /*@
1219*6f263ca3SPeter Brune    SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search.
1220*6f263ca3SPeter Brune 
1221*6f263ca3SPeter Brune    Input Parameters:
1222*6f263ca3SPeter Brune +  linesearch  - linesearch context.
1223*6f263ca3SPeter Brune -  flg  - indicates whether or not to compute norms.
1224*6f263ca3SPeter Brune 
1225*6f263ca3SPeter Brune    Options Database Keys:
1226*6f263ca3SPeter Brune .   -snes_linesearch_norms - turn norm computation on or off.
1227*6f263ca3SPeter Brune 
1228*6f263ca3SPeter Brune    Notes:
1229*6f263ca3SPeter Brune    This is most relevant to the SNES_LINESEARCH_BASIC line search type.
1230*6f263ca3SPeter Brune 
1231*6f263ca3SPeter Brune    Level: intermediate
1232*6f263ca3SPeter Brune 
1233*6f263ca3SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetNorms(), SNESLineSearchComputeNorms(), SNES_LINESEARCH_BASIC
1234*6f263ca3SPeter Brune @*/
1235*6f263ca3SPeter Brune PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg)
1236*6f263ca3SPeter Brune {
1237*6f263ca3SPeter Brune   PetscFunctionBegin;
1238*6f263ca3SPeter Brune   linesearch->norms = flg;
1239*6f263ca3SPeter Brune   PetscFunctionReturn(0);
1240*6f263ca3SPeter Brune }
1241*6f263ca3SPeter Brune 
12426a388c36SPeter Brune #undef __FUNCT__
1243f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetVecs"
1244f40b411bSPeter Brune /*@
1245f1c6b773SPeter Brune    SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context
1246f40b411bSPeter Brune 
1247f40b411bSPeter Brune    Input Parameters:
1248f40b411bSPeter Brune .  linesearch - linesearch context.
1249f40b411bSPeter Brune 
1250f40b411bSPeter Brune    Output Parameters:
1251f40b411bSPeter Brune +  X - The old solution
1252f40b411bSPeter Brune .  F - The old function
1253f40b411bSPeter Brune .  Y - The search direction
1254f40b411bSPeter Brune .  W - The new solution
1255f40b411bSPeter Brune -  G - The new function
1256f40b411bSPeter Brune 
1257f40b411bSPeter Brune    Level: intermediate
1258f40b411bSPeter Brune 
1259f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1260f40b411bSPeter Brune @*/
1261f1c6b773SPeter Brune PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G) {
12626a388c36SPeter Brune   PetscFunctionBegin;
1263f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
12646a388c36SPeter Brune   if (X) {
12656a388c36SPeter Brune     PetscValidPointer(X, 2);
12666a388c36SPeter Brune     *X = linesearch->vec_sol;
12676a388c36SPeter Brune   }
12686a388c36SPeter Brune   if (F) {
12696a388c36SPeter Brune     PetscValidPointer(F, 3);
12706a388c36SPeter Brune     *F = linesearch->vec_func;
12716a388c36SPeter Brune   }
12726a388c36SPeter Brune   if (Y) {
12736a388c36SPeter Brune     PetscValidPointer(Y, 4);
12746a388c36SPeter Brune     *Y = linesearch->vec_update;
12756a388c36SPeter Brune   }
12766a388c36SPeter Brune   if (W) {
12776a388c36SPeter Brune     PetscValidPointer(W, 5);
12786a388c36SPeter Brune     *W = linesearch->vec_sol_new;
12796a388c36SPeter Brune   }
12806a388c36SPeter Brune   if (G) {
12816a388c36SPeter Brune     PetscValidPointer(G, 6);
12826a388c36SPeter Brune     *G = linesearch->vec_func_new;
12836a388c36SPeter Brune   }
12846a388c36SPeter Brune 
12856a388c36SPeter Brune   PetscFunctionReturn(0);
12866a388c36SPeter Brune }
12876a388c36SPeter Brune 
12886a388c36SPeter Brune #undef __FUNCT__
1289f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetVecs"
1290f40b411bSPeter Brune /*@
1291f1c6b773SPeter Brune    SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context
1292f40b411bSPeter Brune 
1293f40b411bSPeter Brune    Input Parameters:
1294f40b411bSPeter Brune +  linesearch - linesearch context.
1295f40b411bSPeter Brune .  X - The old solution
1296f40b411bSPeter Brune .  F - The old function
1297f40b411bSPeter Brune .  Y - The search direction
1298f40b411bSPeter Brune .  W - The new solution
1299f40b411bSPeter Brune -  G - The new function
1300f40b411bSPeter Brune 
1301f40b411bSPeter Brune    Level: intermediate
1302f40b411bSPeter Brune 
1303f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms(), SNESLineSearchGetVecs()
1304f40b411bSPeter Brune @*/
1305f1c6b773SPeter Brune PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G) {
13066a388c36SPeter Brune   PetscFunctionBegin;
1307f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
13086a388c36SPeter Brune   if (X) {
13096a388c36SPeter Brune     PetscValidHeaderSpecific(X,VEC_CLASSID,2);
13106a388c36SPeter Brune     linesearch->vec_sol = X;
13116a388c36SPeter Brune   }
13126a388c36SPeter Brune   if (F) {
13136a388c36SPeter Brune     PetscValidHeaderSpecific(F,VEC_CLASSID,3);
13146a388c36SPeter Brune     linesearch->vec_func = F;
13156a388c36SPeter Brune   }
13166a388c36SPeter Brune   if (Y) {
13176a388c36SPeter Brune     PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
13186a388c36SPeter Brune     linesearch->vec_update = Y;
13196a388c36SPeter Brune   }
13206a388c36SPeter Brune   if (W) {
13216a388c36SPeter Brune     PetscValidHeaderSpecific(W,VEC_CLASSID,5);
13226a388c36SPeter Brune     linesearch->vec_sol_new = W;
13236a388c36SPeter Brune   }
13246a388c36SPeter Brune   if (G) {
13256a388c36SPeter Brune     PetscValidHeaderSpecific(G,VEC_CLASSID,6);
13266a388c36SPeter Brune     linesearch->vec_func_new = G;
13276a388c36SPeter Brune   }
13286a388c36SPeter Brune 
13296a388c36SPeter Brune   PetscFunctionReturn(0);
13306a388c36SPeter Brune }
13316a388c36SPeter Brune 
13326a388c36SPeter Brune #undef __FUNCT__
1333f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchAppendOptionsPrefix"
1334e7058c64SPeter Brune /*@C
1335f1c6b773SPeter Brune    SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1336e7058c64SPeter Brune    SNES options in the database.
1337e7058c64SPeter Brune 
1338cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
1339e7058c64SPeter Brune 
1340e7058c64SPeter Brune    Input Parameters:
1341e7058c64SPeter Brune +  snes - the SNES context
1342e7058c64SPeter Brune -  prefix - the prefix to prepend to all option names
1343e7058c64SPeter Brune 
1344e7058c64SPeter Brune    Notes:
1345e7058c64SPeter Brune    A hyphen (-) must NOT be given at the beginning of the prefix name.
1346e7058c64SPeter Brune    The first character of all runtime options is AUTOMATICALLY the hyphen.
1347e7058c64SPeter Brune 
1348e7058c64SPeter Brune    Level: advanced
1349e7058c64SPeter Brune 
1350f1c6b773SPeter Brune .keywords: SNESLineSearch, append, options, prefix, database
1351e7058c64SPeter Brune 
1352e7058c64SPeter Brune .seealso: SNESGetOptionsPrefix()
1353e7058c64SPeter Brune @*/
1354f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[])
1355e7058c64SPeter Brune {
1356e7058c64SPeter Brune   PetscErrorCode ierr;
1357e7058c64SPeter Brune 
1358e7058c64SPeter Brune   PetscFunctionBegin;
1359f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1360e7058c64SPeter Brune   ierr = PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1361e7058c64SPeter Brune   PetscFunctionReturn(0);
1362e7058c64SPeter Brune }
1363e7058c64SPeter Brune 
1364e7058c64SPeter Brune #undef __FUNCT__
1365f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetOptionsPrefix"
1366e7058c64SPeter Brune /*@C
1367f1c6b773SPeter Brune    SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all
1368f1c6b773SPeter Brune    SNESLineSearch options in the database.
1369e7058c64SPeter Brune 
1370e7058c64SPeter Brune    Not Collective
1371e7058c64SPeter Brune 
1372e7058c64SPeter Brune    Input Parameter:
1373cd7522eaSPeter Brune .  linesearch - the SNESLineSearch context
1374e7058c64SPeter Brune 
1375e7058c64SPeter Brune    Output Parameter:
1376e7058c64SPeter Brune .  prefix - pointer to the prefix string used
1377e7058c64SPeter Brune 
1378e7058c64SPeter Brune    Notes: On the fortran side, the user should pass in a string 'prefix' of
1379e7058c64SPeter Brune    sufficient length to hold the prefix.
1380e7058c64SPeter Brune 
1381e7058c64SPeter Brune    Level: advanced
1382e7058c64SPeter Brune 
1383f1c6b773SPeter Brune .keywords: SNESLineSearch, get, options, prefix, database
1384e7058c64SPeter Brune 
1385e7058c64SPeter Brune .seealso: SNESAppendOptionsPrefix()
1386e7058c64SPeter Brune @*/
1387f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[])
1388e7058c64SPeter Brune {
1389e7058c64SPeter Brune   PetscErrorCode ierr;
1390e7058c64SPeter Brune 
1391e7058c64SPeter Brune   PetscFunctionBegin;
1392f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1393e7058c64SPeter Brune   ierr = PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1394e7058c64SPeter Brune   PetscFunctionReturn(0);
1395e7058c64SPeter Brune }
1396bf7f4e0aSPeter Brune 
1397bf7f4e0aSPeter Brune #undef __FUNCT__
1398f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetWork"
1399f40b411bSPeter Brune /*@
1400f1c6b773SPeter Brune    SNESLineSearchGetWork - Gets work vectors for the line search.
1401f40b411bSPeter Brune 
1402f40b411bSPeter Brune    Input Parameter:
1403f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
1404f40b411bSPeter Brune -  nwork - the number of work vectors
1405f40b411bSPeter Brune 
1406f40b411bSPeter Brune    Level: developer
1407f40b411bSPeter Brune 
1408cd7522eaSPeter Brune    Notes:
1409cd7522eaSPeter Brune    This is typically called at the beginning of a SNESLineSearch or SNESLineSearchShell implementation.
1410cd7522eaSPeter Brune 
1411f1c6b773SPeter Brune .keywords: SNESLineSearch, work, vector
1412f40b411bSPeter Brune 
1413f40b411bSPeter Brune .seealso: SNESDefaultGetWork()
1414f40b411bSPeter Brune @*/
1415f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetWork(SNESLineSearch linesearch, PetscInt nwork)
1416bf7f4e0aSPeter Brune {
1417bf7f4e0aSPeter Brune   PetscErrorCode ierr;
1418bf7f4e0aSPeter Brune   PetscFunctionBegin;
1419bf7f4e0aSPeter Brune   if (linesearch->vec_sol) {
1420bf7f4e0aSPeter Brune     ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr);
1421bf7f4e0aSPeter Brune   } else {
1422bf7f4e0aSPeter Brune     SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
1423bf7f4e0aSPeter Brune   }
1424bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1425bf7f4e0aSPeter Brune }
1426bf7f4e0aSPeter Brune 
14276a388c36SPeter Brune 
1428bf7f4e0aSPeter Brune #undef __FUNCT__
1429f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetSuccess"
1430f40b411bSPeter Brune /*@
1431f1c6b773SPeter Brune    SNESLineSearchGetSuccess - Gets the success/failure status of the last line search application
1432f40b411bSPeter Brune 
1433f40b411bSPeter Brune    Input Parameters:
1434f40b411bSPeter Brune .  linesearch - linesearch context.
1435f40b411bSPeter Brune 
1436f40b411bSPeter Brune    Output Parameters:
1437f40b411bSPeter Brune .  success - The success or failure status.
1438f40b411bSPeter Brune 
1439cd7522eaSPeter Brune    Notes:
1440cd7522eaSPeter Brune    This is typically called after SNESLineSearchApply in order to determine if the line-search failed
1441cd7522eaSPeter Brune    (and set the SNES convergence accordingly).
1442cd7522eaSPeter Brune 
1443f40b411bSPeter Brune    Level: intermediate
1444f40b411bSPeter Brune 
1445f1c6b773SPeter Brune .seealso: SNESLineSearchSetSuccess()
1446f40b411bSPeter Brune @*/
1447f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetSuccess(SNESLineSearch linesearch, PetscBool *success)
1448bf7f4e0aSPeter Brune {
1449bf7f4e0aSPeter Brune   PetscFunctionBegin;
1450f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
14516a388c36SPeter Brune   PetscValidPointer(success, 2);
1452bf7f4e0aSPeter Brune   if (success) {
1453bf7f4e0aSPeter Brune     *success = linesearch->success;
1454bf7f4e0aSPeter Brune   }
1455bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1456bf7f4e0aSPeter Brune }
1457bf7f4e0aSPeter Brune 
1458bf7f4e0aSPeter Brune #undef __FUNCT__
1459f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetSuccess"
1460f40b411bSPeter Brune /*@
1461f1c6b773SPeter Brune    SNESLineSearchSetSuccess - Sets the success/failure status of the last line search application
1462f40b411bSPeter Brune 
1463f40b411bSPeter Brune    Input Parameters:
1464f40b411bSPeter Brune +  linesearch - linesearch context.
1465f40b411bSPeter Brune -  success - The success or failure status.
1466f40b411bSPeter Brune 
1467cd7522eaSPeter Brune    Notes:
1468cd7522eaSPeter Brune    This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set
1469cd7522eaSPeter Brune    the success or failure of the line search method.
1470cd7522eaSPeter Brune 
1471f40b411bSPeter Brune    Level: intermediate
1472f40b411bSPeter Brune 
1473f1c6b773SPeter Brune .seealso: SNESLineSearchGetSuccess()
1474f40b411bSPeter Brune @*/
1475f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetSuccess(SNESLineSearch linesearch, PetscBool success)
14766a388c36SPeter Brune {
1477f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
14786a388c36SPeter Brune   PetscFunctionBegin;
14796a388c36SPeter Brune   linesearch->success = success;
14806a388c36SPeter Brune   PetscFunctionReturn(0);
14816a388c36SPeter Brune }
14826a388c36SPeter Brune 
14836a388c36SPeter Brune #undef __FUNCT__
1484f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetVIFunctions"
14859bd66eb0SPeter Brune /*@C
1486f1c6b773SPeter Brune    SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.
14879bd66eb0SPeter Brune 
14889bd66eb0SPeter Brune    Input Parameters:
14899bd66eb0SPeter Brune +  snes - nonlinear context obtained from SNESCreate()
14909bd66eb0SPeter Brune .  projectfunc - function for projecting the function to the bounds
14919bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
14929bd66eb0SPeter Brune 
14939bd66eb0SPeter Brune    Logically Collective on SNES
14949bd66eb0SPeter Brune 
14959bd66eb0SPeter Brune    Calling sequence of projectfunc:
14969bd66eb0SPeter Brune .vb
14979bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X)
14989bd66eb0SPeter Brune .ve
14999bd66eb0SPeter Brune 
15009bd66eb0SPeter Brune     Input parameters for projectfunc:
15019bd66eb0SPeter Brune +   snes - nonlinear context
15029bd66eb0SPeter Brune -   X - current solution
15039bd66eb0SPeter Brune 
1504cd7522eaSPeter Brune     Output parameters for projectfunc:
15059bd66eb0SPeter Brune .   X - Projected solution
15069bd66eb0SPeter Brune 
15079bd66eb0SPeter Brune    Calling sequence of normfunc:
15089bd66eb0SPeter Brune .vb
15099bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm)
15109bd66eb0SPeter Brune .ve
15119bd66eb0SPeter Brune 
1512cd7522eaSPeter Brune     Input parameters for normfunc:
15139bd66eb0SPeter Brune +   snes - nonlinear context
15149bd66eb0SPeter Brune .   X - current solution
15159bd66eb0SPeter Brune -   F - current residual
15169bd66eb0SPeter Brune 
1517cd7522eaSPeter Brune     Output parameters for normfunc:
15189bd66eb0SPeter Brune .   fnorm - VI-specific norm of the function
15199bd66eb0SPeter Brune 
1520cd7522eaSPeter Brune     Notes:
1521cd7522eaSPeter Brune     The VI solvers require projection of the solution to the feasible set.  projectfunc should implement this.
1522cd7522eaSPeter Brune 
1523cd7522eaSPeter Brune     The VI solvers require special evaluation of the function norm such that the norm is only calculated
1524cd7522eaSPeter Brune     on the inactive set.  This should be implemented by normfunc.
15259bd66eb0SPeter Brune 
15269bd66eb0SPeter Brune     Level: developer
15279bd66eb0SPeter Brune 
15289bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, set, line search
15299bd66eb0SPeter Brune 
1530f1c6b773SPeter Brune .seealso: SNESLineSearchGetVIFunctions(), SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck()
15319bd66eb0SPeter Brune @*/
1532f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc)
15339bd66eb0SPeter Brune {
15349bd66eb0SPeter Brune   PetscFunctionBegin;
1535f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
15369bd66eb0SPeter Brune   if (projectfunc) linesearch->ops->viproject = projectfunc;
15379bd66eb0SPeter Brune   if (normfunc) linesearch->ops->vinorm = normfunc;
15389bd66eb0SPeter Brune   PetscFunctionReturn(0);
15399bd66eb0SPeter Brune }
15409bd66eb0SPeter Brune 
15419bd66eb0SPeter Brune #undef __FUNCT__
1542f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetVIFunctions"
15439bd66eb0SPeter Brune /*@C
1544f1c6b773SPeter Brune    SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.
15459bd66eb0SPeter Brune 
15469bd66eb0SPeter Brune    Input Parameters:
15479bd66eb0SPeter Brune .  snes - nonlinear context obtained from SNESCreate()
15489bd66eb0SPeter Brune 
15499bd66eb0SPeter Brune    Output Parameters:
15509bd66eb0SPeter Brune +  projectfunc - function for projecting the function to the bounds
15519bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
15529bd66eb0SPeter Brune 
15539bd66eb0SPeter Brune    Logically Collective on SNES
15549bd66eb0SPeter Brune 
15559bd66eb0SPeter Brune     Level: developer
15569bd66eb0SPeter Brune 
15579bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, get, line search
15589bd66eb0SPeter Brune 
1559f1c6b773SPeter Brune .seealso: SNESLineSearchSetVIFunctions(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
15609bd66eb0SPeter Brune @*/
1561f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc)
15629bd66eb0SPeter Brune {
15639bd66eb0SPeter Brune   PetscFunctionBegin;
15649bd66eb0SPeter Brune   if (projectfunc) *projectfunc = linesearch->ops->viproject;
15659bd66eb0SPeter Brune   if (normfunc) *normfunc = linesearch->ops->vinorm;
15669bd66eb0SPeter Brune   PetscFunctionReturn(0);
15679bd66eb0SPeter Brune }
15689bd66eb0SPeter Brune 
15699bd66eb0SPeter Brune #undef __FUNCT__
1570f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchRegister"
1571bf7f4e0aSPeter Brune /*@C
1572f1c6b773SPeter Brune   SNESLineSearchRegister - See SNESLineSearchRegisterDynamic()
1573bf7f4e0aSPeter Brune 
1574bf7f4e0aSPeter Brune   Level: advanced
1575bf7f4e0aSPeter Brune @*/
1576f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNESLineSearch))
1577bf7f4e0aSPeter Brune {
1578bf7f4e0aSPeter Brune   char           fullname[PETSC_MAX_PATH_LEN];
1579bf7f4e0aSPeter Brune   PetscErrorCode ierr;
1580bf7f4e0aSPeter Brune 
1581bf7f4e0aSPeter Brune   PetscFunctionBegin;
1582bf7f4e0aSPeter Brune   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
1583f1c6b773SPeter Brune   ierr = PetscFListAdd(&SNESLineSearchList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
1584bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1585bf7f4e0aSPeter Brune }
1586