xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision aaa7dc30da3270cff6cb10b1db605b2ca746f216)
1b45d2f2cSJed Brown #include <petsc-private/linesearchimpl.h> /*I "petscsnes.h" I*/
2bf7f4e0aSPeter Brune 
3f1c6b773SPeter Brune PetscBool         SNESLineSearchRegisterAllCalled = PETSC_FALSE;
40298fd71SBarry Smith PetscFunctionList SNESLineSearchList              = 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 
29bf388a1fSBarry Smith PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch)
30bf388a1fSBarry Smith {
31bf7f4e0aSPeter Brune   PetscErrorCode ierr;
32f1c6b773SPeter Brune   SNESLineSearch linesearch;
33bf388a1fSBarry Smith 
34bf7f4e0aSPeter Brune   PetscFunctionBegin;
35ea5d4fccSPeter Brune   PetscValidPointer(outlinesearch,2);
36f34a81feSMatthew G. Knepley   ierr = SNESInitializePackage();CHKERRQ(ierr);
370298fd71SBarry Smith   *outlinesearch = NULL;
38f5af7f23SKarl Rupp 
3967c2884eSBarry Smith   ierr = PetscHeaderCreate(linesearch,_p_LineSearch,struct _LineSearchOps,SNESLINESEARCH_CLASSID, "SNESLineSearch","Linesearch","SNESLineSearch",comm,SNESLineSearchDestroy,SNESLineSearchView);CHKERRQ(ierr);
40bf7f4e0aSPeter Brune 
410298fd71SBarry Smith   linesearch->ops->precheck  = NULL;
420298fd71SBarry Smith   linesearch->ops->postcheck = NULL;
43bf7f4e0aSPeter Brune 
440298fd71SBarry Smith   linesearch->vec_sol_new  = NULL;
450298fd71SBarry Smith   linesearch->vec_func_new = NULL;
460298fd71SBarry Smith   linesearch->vec_sol      = NULL;
470298fd71SBarry Smith   linesearch->vec_func     = NULL;
480298fd71SBarry Smith   linesearch->vec_update   = NULL;
499bd66eb0SPeter Brune 
50bf7f4e0aSPeter Brune   linesearch->lambda       = 1.0;
51bf7f4e0aSPeter Brune   linesearch->fnorm        = 1.0;
52bf7f4e0aSPeter Brune   linesearch->ynorm        = 1.0;
53bf7f4e0aSPeter Brune   linesearch->xnorm        = 1.0;
54bf7f4e0aSPeter Brune   linesearch->success      = PETSC_TRUE;
55bf7f4e0aSPeter Brune   linesearch->norms        = PETSC_TRUE;
56bf7f4e0aSPeter Brune   linesearch->keeplambda   = PETSC_FALSE;
57bf7f4e0aSPeter Brune   linesearch->damping      = 1.0;
58bf7f4e0aSPeter Brune   linesearch->maxstep      = 1e8;
59bf7f4e0aSPeter Brune   linesearch->steptol      = 1e-12;
60516fe3c3SPeter Brune   linesearch->rtol         = 1e-8;
61516fe3c3SPeter Brune   linesearch->atol         = 1e-15;
62516fe3c3SPeter Brune   linesearch->ltol         = 1e-8;
630298fd71SBarry Smith   linesearch->precheckctx  = NULL;
640298fd71SBarry Smith   linesearch->postcheckctx = NULL;
6559405d5eSPeter Brune   linesearch->max_its      = 1;
66bf7f4e0aSPeter Brune   linesearch->setupcalled  = PETSC_FALSE;
67bf7f4e0aSPeter Brune   *outlinesearch           = linesearch;
68bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
69bf7f4e0aSPeter Brune }
70bf7f4e0aSPeter Brune 
71bf7f4e0aSPeter Brune #undef __FUNCT__
72f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetUp"
73f40b411bSPeter Brune /*@
7478bcb3b5SPeter Brune    SNESLineSearchSetUp - Prepares the line search for being applied by allocating
7578bcb3b5SPeter Brune    any required vectors.
76f40b411bSPeter Brune 
77cd7522eaSPeter Brune    Collective on SNESLineSearch
78f40b411bSPeter Brune 
79f40b411bSPeter Brune    Input Parameters:
80f40b411bSPeter Brune .  linesearch - The LineSearch instance.
81f40b411bSPeter Brune 
82cd7522eaSPeter Brune    Notes:
83cd7522eaSPeter Brune    For most cases, this needn't be called outside of SNESLineSearchApply().
84cd7522eaSPeter Brune    The only current case where this is called outside of this is for the VI
8578bcb3b5SPeter Brune    solvers, which modify the solution and work vectors before the first call
86cd7522eaSPeter Brune    of SNESLineSearchApply, requiring the SNESLineSearch work vectors to be
87cd7522eaSPeter Brune    allocated upfront.
88cd7522eaSPeter Brune 
8978bcb3b5SPeter Brune    Level: advanced
90f40b411bSPeter Brune 
91f1c6b773SPeter Brune .keywords: SNESLineSearch, SetUp
92f40b411bSPeter Brune 
93f1c6b773SPeter Brune .seealso: SNESLineSearchReset()
94f40b411bSPeter Brune @*/
95f40b411bSPeter Brune 
96bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch)
97bf388a1fSBarry Smith {
98bf7f4e0aSPeter Brune   PetscErrorCode ierr;
99bf388a1fSBarry Smith 
100bf7f4e0aSPeter Brune   PetscFunctionBegin;
101bf7f4e0aSPeter Brune   if (!((PetscObject)linesearch)->type_name) {
1021a4f838cSPeter Brune     ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr);
103bf7f4e0aSPeter Brune   }
104bf7f4e0aSPeter Brune   if (!linesearch->setupcalled) {
1059bd66eb0SPeter Brune     if (!linesearch->vec_sol_new) {
106bf7f4e0aSPeter Brune       ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new);CHKERRQ(ierr);
1079bd66eb0SPeter Brune     }
1089bd66eb0SPeter Brune     if (!linesearch->vec_func_new) {
1099bd66eb0SPeter Brune       ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new);CHKERRQ(ierr);
1109bd66eb0SPeter Brune     }
111bf7f4e0aSPeter Brune     if (linesearch->ops->setup) {
112bf7f4e0aSPeter Brune       ierr = (*linesearch->ops->setup)(linesearch);CHKERRQ(ierr);
113bf7f4e0aSPeter Brune     }
114ed07d7d7SPeter Brune     if (!linesearch->ops->snesfunc) {ierr = SNESLineSearchSetFunction(linesearch,SNESComputeFunction);CHKERRQ(ierr);}
115bf7f4e0aSPeter Brune     linesearch->lambda      = linesearch->damping;
116bf7f4e0aSPeter Brune     linesearch->setupcalled = PETSC_TRUE;
117bf7f4e0aSPeter Brune   }
118bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
119bf7f4e0aSPeter Brune }
120bf7f4e0aSPeter Brune 
121bf7f4e0aSPeter Brune #undef __FUNCT__
122f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchReset"
123f40b411bSPeter Brune 
124f40b411bSPeter Brune /*@
12578bcb3b5SPeter Brune    SNESLineSearchReset - Undoes the SetUp and deletes any Vecs or Mats allocated by the line search.
126f40b411bSPeter Brune 
127f1c6b773SPeter Brune    Collective on SNESLineSearch
128f40b411bSPeter Brune 
129f40b411bSPeter Brune    Input Parameters:
130f40b411bSPeter Brune .  linesearch - The LineSearch instance.
131f40b411bSPeter Brune 
13278bcb3b5SPeter Brune    Level: intermediate
133f40b411bSPeter Brune 
134cd7522eaSPeter Brune .keywords: SNESLineSearch, Reset
135f40b411bSPeter Brune 
136f1c6b773SPeter Brune .seealso: SNESLineSearchSetUp()
137f40b411bSPeter Brune @*/
138f40b411bSPeter Brune 
139bf388a1fSBarry Smith PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch)
140bf388a1fSBarry Smith {
141bf7f4e0aSPeter Brune   PetscErrorCode ierr;
142bf388a1fSBarry Smith 
143bf7f4e0aSPeter Brune   PetscFunctionBegin;
144f5af7f23SKarl Rupp   if (linesearch->ops->reset) (*linesearch->ops->reset)(linesearch);
145f5af7f23SKarl Rupp 
146bf7f4e0aSPeter Brune   ierr = VecDestroy(&linesearch->vec_sol_new);CHKERRQ(ierr);
147bf7f4e0aSPeter Brune   ierr = VecDestroy(&linesearch->vec_func_new);CHKERRQ(ierr);
148bf7f4e0aSPeter Brune 
149bf7f4e0aSPeter Brune   ierr = VecDestroyVecs(linesearch->nwork, &linesearch->work);CHKERRQ(ierr);
150f5af7f23SKarl Rupp 
151bf7f4e0aSPeter Brune   linesearch->nwork       = 0;
152bf7f4e0aSPeter Brune   linesearch->setupcalled = PETSC_FALSE;
153bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
154bf7f4e0aSPeter Brune }
155bf7f4e0aSPeter Brune 
156ed07d7d7SPeter Brune #undef __FUNCT__
157ed07d7d7SPeter Brune #define __FUNCT__ "SNESLineSearchSetFunction"
158ed07d7d7SPeter Brune /*@C
159ed07d7d7SPeter Brune    SNESLineSearchSetPreCheck - Sets the function evaluation context on the
160ed07d7d7SPeter Brune 
161ed07d7d7SPeter Brune    Input Parameters:
162ed07d7d7SPeter Brune .  linesearch - the SNESLineSearch context
163ed07d7d7SPeter Brune 
164ed07d7d7SPeter Brune    Output Parameters:
165ed07d7d7SPeter Brune +  func       - [optional] function evaluation routine
166ed07d7d7SPeter Brune 
167ed07d7d7SPeter Brune    Level: developer
168ed07d7d7SPeter Brune 
169ed07d7d7SPeter Brune .keywords: get, linesearch, pre-check
170ed07d7d7SPeter Brune 
171ed07d7d7SPeter Brune .seealso: SNESLineSearchGetPostCheck(), SNESLineSearchSetPreCheck()
172ed07d7d7SPeter Brune @*/
173ed07d7d7SPeter Brune PetscErrorCode  SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES,Vec,Vec))
174ed07d7d7SPeter Brune {
175ed07d7d7SPeter Brune   PetscFunctionBegin;
176ed07d7d7SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
177ed07d7d7SPeter Brune   linesearch->ops->snesfunc = func;
178ed07d7d7SPeter Brune   PetscFunctionReturn(0);
179ed07d7d7SPeter Brune }
180ed07d7d7SPeter Brune 
181ed07d7d7SPeter Brune 
1826b2b7091SBarry Smith /*MC
1836b2b7091SBarry Smith     SNESLineSearchPreCheckFunction - functional form passed to check before line search is called
1846b2b7091SBarry Smith 
1856b2b7091SBarry Smith      Synopsis:
186*aaa7dc30SBarry Smith      #include <petscsnes.h>
1876b2b7091SBarry Smith      SNESLineSearchPreCheckFunction(SNESLineSearch snes,Vec x,Vec y, PetscBool *changed);
1886b2b7091SBarry Smith 
1896b2b7091SBarry Smith        Input Parameters:
1906b2b7091SBarry Smith +      x - solution vector
1916b2b7091SBarry Smith .      y - search direction vector
1926b2b7091SBarry Smith -      changed - flag to indicate the precheck changed x or y.
1936b2b7091SBarry Smith 
194878cb397SSatish Balay    Level: advanced
195878cb397SSatish Balay 
1966b2b7091SBarry Smith .seealso:   SNESLineSearchSetPreCheck(), SNESLineSearchSetPostCheck()
1976b2b7091SBarry Smith M*/
1986b2b7091SBarry Smith 
19986d74e61SPeter Brune #undef __FUNCT__
200f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetPreCheck"
20186d74e61SPeter Brune /*@C
202f1c6b773SPeter Brune    SNESLineSearchSetPreCheck - Sets a pre-check function for the line search routine.
20386d74e61SPeter Brune 
204f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
20586d74e61SPeter Brune 
20686d74e61SPeter Brune    Input Parameters:
207f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
2086b2b7091SBarry Smith .  SNESLineSearchPreCheckFunction - [optional] function evaluation routine
20986d74e61SPeter Brune -  ctx        - [optional] user-defined context for private data for the
2100298fd71SBarry Smith                 function evaluation routine (may be NULL)
21186d74e61SPeter Brune 
21286d74e61SPeter Brune 
21386d74e61SPeter Brune    Level: intermediate
21486d74e61SPeter Brune 
21586d74e61SPeter Brune .keywords: set, linesearch, pre-check
21686d74e61SPeter Brune 
217f1c6b773SPeter Brune .seealso: SNESLineSearchSetPostCheck()
21886d74e61SPeter Brune @*/
2196b2b7091SBarry Smith PetscErrorCode  SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*SNESLineSearchPreCheckFunction)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void *ctx)
22086d74e61SPeter Brune {
2219bd66eb0SPeter Brune   PetscFunctionBegin;
222f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
2236b2b7091SBarry Smith   if (SNESLineSearchPreCheckFunction) linesearch->ops->precheck = SNESLineSearchPreCheckFunction;
22486d74e61SPeter Brune   if (ctx) linesearch->precheckctx = ctx;
22586d74e61SPeter Brune   PetscFunctionReturn(0);
22686d74e61SPeter Brune }
22786d74e61SPeter Brune 
22886d74e61SPeter Brune #undef __FUNCT__
229f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetPreCheck"
23086d74e61SPeter Brune /*@C
231cd7522eaSPeter Brune    SNESLineSearchSetPreCheck - Gets the pre-check function for the line search routine.
23286d74e61SPeter Brune 
23386d74e61SPeter Brune    Input Parameters:
234f1c6b773SPeter Brune .  linesearch - the SNESLineSearch context
23586d74e61SPeter Brune 
23686d74e61SPeter Brune    Output Parameters:
23786d74e61SPeter Brune +  func       - [optional] function evaluation routine
23886d74e61SPeter Brune -  ctx        - [optional] user-defined context for private data for the
2390298fd71SBarry Smith                 function evaluation routine (may be NULL)
24086d74e61SPeter Brune 
24186d74e61SPeter Brune    Level: intermediate
24286d74e61SPeter Brune 
24386d74e61SPeter Brune .keywords: get, linesearch, pre-check
24486d74e61SPeter Brune 
245f1c6b773SPeter Brune .seealso: SNESLineSearchGetPostCheck(), SNESLineSearchSetPreCheck()
24686d74e61SPeter Brune @*/
2476b2b7091SBarry Smith PetscErrorCode  SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**SNESLineSearchPreCheckFunction)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void **ctx)
24886d74e61SPeter Brune {
2499bd66eb0SPeter Brune   PetscFunctionBegin;
250f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
2516b2b7091SBarry Smith   if (SNESLineSearchPreCheckFunction) *SNESLineSearchPreCheckFunction = linesearch->ops->precheck;
25286d74e61SPeter Brune   if (ctx) *ctx = linesearch->precheckctx;
25386d74e61SPeter Brune   PetscFunctionReturn(0);
25486d74e61SPeter Brune }
25586d74e61SPeter Brune 
2566b2b7091SBarry Smith /*MC
2576b2b7091SBarry Smith     SNESLineSearchPostheckFunction - functional form that is called after line search is complete
2586b2b7091SBarry Smith 
2596b2b7091SBarry Smith      Synopsis:
260*aaa7dc30SBarry Smith      #include <petscsnes.h>
2616b2b7091SBarry Smith      SNESLineSearchPostheckFunction(SNESLineSearch linesearch,Vec x,Vec y,  Vec w, *changed_y, PetscBool *changed_w);
2626b2b7091SBarry Smith 
2636b2b7091SBarry Smith      Input Parameters:
2646b2b7091SBarry Smith +      x - old solution vector
2656b2b7091SBarry Smith .      y - search direction vector
2666b2b7091SBarry Smith .      w - new solution vector
2676b2b7091SBarry Smith .      changed_y - indicates that the line search changed y
2686b2b7091SBarry Smith -      changed_w - indicates that the line search changed w
2696b2b7091SBarry Smith 
270878cb397SSatish Balay    Level: advanced
2716b2b7091SBarry Smith 
2726b2b7091SBarry Smith .seealso:   SNESLineSearchSetPreCheck(), SNESLineSearchSetPostCheck()
2736b2b7091SBarry Smith M*/
2746b2b7091SBarry Smith 
27586d74e61SPeter Brune #undef __FUNCT__
276f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetPostCheck"
27786d74e61SPeter Brune /*@C
278f1c6b773SPeter Brune    SNESLineSearchSetPostCheck - Sets a post-check function for the line search routine.
27986d74e61SPeter Brune 
280f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
28186d74e61SPeter Brune 
28286d74e61SPeter Brune    Input Parameters:
283f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
2846b2b7091SBarry Smith .  SNESLineSearchPostCheckFunction - [optional] function evaluation routine
28586d74e61SPeter Brune -  ctx        - [optional] user-defined context for private data for the
2860298fd71SBarry Smith                 function evaluation routine (may be NULL)
28786d74e61SPeter Brune 
28886d74e61SPeter Brune    Level: intermediate
28986d74e61SPeter Brune 
29086d74e61SPeter Brune .keywords: set, linesearch, post-check
29186d74e61SPeter Brune 
292f1c6b773SPeter Brune .seealso: SNESLineSearchSetPreCheck()
29386d74e61SPeter Brune @*/
2946b2b7091SBarry Smith PetscErrorCode  SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*SNESLineSearchPostCheckFunction)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void *ctx)
29586d74e61SPeter Brune {
29686d74e61SPeter Brune   PetscFunctionBegin;
297f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
2986b2b7091SBarry Smith   if (SNESLineSearchPostCheckFunction) linesearch->ops->postcheck = SNESLineSearchPostCheckFunction;
29986d74e61SPeter Brune   if (ctx) linesearch->postcheckctx = ctx;
30086d74e61SPeter Brune   PetscFunctionReturn(0);
30186d74e61SPeter Brune }
30286d74e61SPeter Brune 
30386d74e61SPeter Brune #undef __FUNCT__
304f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetPostCheck"
30586d74e61SPeter Brune /*@C
306f1c6b773SPeter Brune    SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine.
30786d74e61SPeter Brune 
30886d74e61SPeter Brune    Input Parameters:
309f1c6b773SPeter Brune .  linesearch - the SNESLineSearch context
31086d74e61SPeter Brune 
31186d74e61SPeter Brune    Output Parameters:
3126b2b7091SBarry Smith +  SNESLineSearchPostCheckFunction - [optional] function evaluation routine
31386d74e61SPeter Brune -  ctx        - [optional] user-defined context for private data for the
3140298fd71SBarry Smith                 function evaluation routine (may be NULL)
31586d74e61SPeter Brune 
31686d74e61SPeter Brune    Level: intermediate
31786d74e61SPeter Brune 
31886d74e61SPeter Brune .keywords: get, linesearch, post-check
31986d74e61SPeter Brune 
320f1c6b773SPeter Brune .seealso: SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck()
32186d74e61SPeter Brune @*/
3226b2b7091SBarry Smith PetscErrorCode  SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**SNESLineSearchPostCheckFunction)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void **ctx)
32386d74e61SPeter Brune {
3249bd66eb0SPeter Brune   PetscFunctionBegin;
325f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
3266b2b7091SBarry Smith   if (SNESLineSearchPostCheckFunction) *SNESLineSearchPostCheckFunction = linesearch->ops->postcheck;
32786d74e61SPeter Brune   if (ctx) *ctx = linesearch->postcheckctx;
32886d74e61SPeter Brune   PetscFunctionReturn(0);
32986d74e61SPeter Brune }
33086d74e61SPeter Brune 
331bf7f4e0aSPeter Brune #undef __FUNCT__
332f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPreCheck"
333f40b411bSPeter Brune /*@
334f1c6b773SPeter Brune    SNESLineSearchPreCheck - Prepares the line search for being applied.
335f40b411bSPeter Brune 
336cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
337f40b411bSPeter Brune 
338f40b411bSPeter Brune    Input Parameters:
3397b1df9c1SPeter Brune +  linesearch - The linesearch instance.
3407b1df9c1SPeter Brune .  X - The current solution
3417b1df9c1SPeter Brune -  Y - The step direction
342f40b411bSPeter Brune 
343f40b411bSPeter Brune    Output Parameters:
3448e557f58SPeter Brune .  changed - Indicator that the precheck routine has changed anything
345f40b411bSPeter Brune 
346f40b411bSPeter Brune    Level: Beginner
347f40b411bSPeter Brune 
348f1c6b773SPeter Brune .keywords: SNESLineSearch, Create
349f40b411bSPeter Brune 
350f1c6b773SPeter Brune .seealso: SNESLineSearchPostCheck()
351f40b411bSPeter Brune @*/
3527b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed)
353bf7f4e0aSPeter Brune {
354bf7f4e0aSPeter Brune   PetscErrorCode ierr;
3555fd66863SKarl Rupp 
356bf7f4e0aSPeter Brune   PetscFunctionBegin;
357bf7f4e0aSPeter Brune   *changed = PETSC_FALSE;
3586b2b7091SBarry Smith   if (linesearch->ops->precheck) {
3596b2b7091SBarry Smith     ierr = (*linesearch->ops->precheck)(linesearch, X, Y, changed, linesearch->precheckctx);CHKERRQ(ierr);
36038bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed,4);
361bf7f4e0aSPeter Brune   }
362bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
363bf7f4e0aSPeter Brune }
364bf7f4e0aSPeter Brune 
365bf7f4e0aSPeter Brune #undef __FUNCT__
366f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPostCheck"
367f40b411bSPeter Brune /*@
368f1c6b773SPeter Brune    SNESLineSearchPostCheck - Prepares the line search for being applied.
369f40b411bSPeter Brune 
370cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
371f40b411bSPeter Brune 
372f40b411bSPeter Brune    Input Parameters:
3737b1df9c1SPeter Brune +  linesearch - The linesearch context
3747b1df9c1SPeter Brune .  X - The last solution
3757b1df9c1SPeter Brune .  Y - The step direction
3767b1df9c1SPeter Brune -  W - The updated solution, W = X + lambda*Y for some lambda
377f40b411bSPeter Brune 
378f40b411bSPeter Brune    Output Parameters:
37978bcb3b5SPeter Brune +  changed_Y - Indicator if the direction Y has been changed.
38078bcb3b5SPeter Brune -  changed_W - Indicator if the new candidate solution W has been changed.
381f40b411bSPeter Brune 
382f40b411bSPeter Brune    Level: Intermediate
383f40b411bSPeter Brune 
384f1c6b773SPeter Brune .keywords: SNESLineSearch, Create
385f40b411bSPeter Brune 
386f1c6b773SPeter Brune .seealso: SNESLineSearchPreCheck()
387f40b411bSPeter Brune @*/
3887b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch,Vec X,Vec Y,Vec W,PetscBool *changed_Y,PetscBool *changed_W)
389bf7f4e0aSPeter Brune {
390bf7f4e0aSPeter Brune   PetscErrorCode ierr;
391bf388a1fSBarry Smith 
392bf7f4e0aSPeter Brune   PetscFunctionBegin;
393bf7f4e0aSPeter Brune   *changed_Y = PETSC_FALSE;
394bf7f4e0aSPeter Brune   *changed_W = PETSC_FALSE;
3956b2b7091SBarry Smith   if (linesearch->ops->postcheck) {
3966b2b7091SBarry Smith     ierr = (*linesearch->ops->postcheck)(linesearch,X,Y,W,changed_Y,changed_W,linesearch->postcheckctx);CHKERRQ(ierr);
39738bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed_Y,5);
39838bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed_W,6);
39986d74e61SPeter Brune   }
40086d74e61SPeter Brune   PetscFunctionReturn(0);
40186d74e61SPeter Brune }
40286d74e61SPeter Brune 
40386d74e61SPeter Brune #undef __FUNCT__
404f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPreCheckPicard"
40586d74e61SPeter Brune /*@C
40686d74e61SPeter Brune    SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration
40786d74e61SPeter Brune 
408cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
40986d74e61SPeter Brune 
41086d74e61SPeter Brune    Input Arguments:
41186d74e61SPeter Brune +  linesearch - linesearch context
41286d74e61SPeter Brune .  X - base state for this step
41386d74e61SPeter Brune .  Y - initial correction
41486d74e61SPeter Brune 
41586d74e61SPeter Brune    Output Arguments:
41686d74e61SPeter Brune +  Y - correction, possibly modified
41786d74e61SPeter Brune -  changed - flag indicating that Y was modified
41886d74e61SPeter Brune 
41986d74e61SPeter Brune    Options Database Key:
420cd7522eaSPeter Brune +  -snes_linesearch_precheck_picard - activate this routine
421cd7522eaSPeter Brune -  -snes_linesearch_precheck_picard_angle - angle
42286d74e61SPeter Brune 
42386d74e61SPeter Brune    Level: advanced
42486d74e61SPeter Brune 
42586d74e61SPeter Brune    Notes:
42686d74e61SPeter Brune    This function should be passed to SNESLineSearchSetPreCheck()
42786d74e61SPeter Brune 
42886d74e61SPeter Brune    The justification for this method involves the linear convergence of a Picard iteration
42986d74e61SPeter Brune    so the Picard linearization should be provided in place of the "Jacobian". This correction
43086d74e61SPeter Brune    is generally not useful when using a Newton linearization.
43186d74e61SPeter Brune 
43286d74e61SPeter Brune    Reference:
43386d74e61SPeter Brune    Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology.
43486d74e61SPeter Brune 
43586d74e61SPeter Brune .seealso: SNESLineSearchSetPreCheck()
43686d74e61SPeter Brune @*/
437f1c6b773SPeter Brune PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed,void *ctx)
43886d74e61SPeter Brune {
43986d74e61SPeter Brune   PetscErrorCode ierr;
44086d74e61SPeter Brune   PetscReal      angle = *(PetscReal*)linesearch->precheckctx;
44186d74e61SPeter Brune   Vec            Ylast;
44286d74e61SPeter Brune   PetscScalar    dot;
44386d74e61SPeter Brune   PetscInt       iter;
44486d74e61SPeter Brune   PetscReal      ynorm,ylastnorm,theta,angle_radians;
44586d74e61SPeter Brune   SNES           snes;
44686d74e61SPeter Brune 
44786d74e61SPeter Brune   PetscFunctionBegin;
448f1c6b773SPeter Brune   ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
44986d74e61SPeter Brune   ierr = PetscObjectQuery((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject*)&Ylast);CHKERRQ(ierr);
45086d74e61SPeter Brune   if (!Ylast) {
45186d74e61SPeter Brune     ierr = VecDuplicate(Y,&Ylast);CHKERRQ(ierr);
45286d74e61SPeter Brune     ierr = PetscObjectCompose((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject)Ylast);CHKERRQ(ierr);
45386d74e61SPeter Brune     ierr = PetscObjectDereference((PetscObject)Ylast);CHKERRQ(ierr);
45486d74e61SPeter Brune   }
45586d74e61SPeter Brune   ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr);
45686d74e61SPeter Brune   if (iter < 2) {
45786d74e61SPeter Brune     ierr     = VecCopy(Y,Ylast);CHKERRQ(ierr);
45886d74e61SPeter Brune     *changed = PETSC_FALSE;
45986d74e61SPeter Brune     PetscFunctionReturn(0);
46086d74e61SPeter Brune   }
46186d74e61SPeter Brune 
46286d74e61SPeter Brune   ierr = VecDot(Y,Ylast,&dot);CHKERRQ(ierr);
46386d74e61SPeter Brune   ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr);
46486d74e61SPeter Brune   ierr = VecNorm(Ylast,NORM_2,&ylastnorm);CHKERRQ(ierr);
46586d74e61SPeter Brune   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
46686d74e61SPeter Brune   theta         = acos((double)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0));
46786d74e61SPeter Brune   angle_radians = angle * PETSC_PI / 180.;
46886d74e61SPeter Brune   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
46986d74e61SPeter Brune     /* Modify the step Y */
47086d74e61SPeter Brune     PetscReal alpha,ydiffnorm;
47186d74e61SPeter Brune     ierr  = VecAXPY(Ylast,-1.0,Y);CHKERRQ(ierr);
47286d74e61SPeter Brune     ierr  = VecNorm(Ylast,NORM_2,&ydiffnorm);CHKERRQ(ierr);
47386d74e61SPeter Brune     alpha = ylastnorm / ydiffnorm;
47486d74e61SPeter Brune     ierr  = VecCopy(Y,Ylast);CHKERRQ(ierr);
47586d74e61SPeter Brune     ierr  = VecScale(Y,alpha);CHKERRQ(ierr);
476c69d1a72SBarry Smith     ierr  = PetscInfo3(snes,"Angle %14.12e degrees less than threshold %14.12e, corrected step by alpha=%14.12e\n",(double)(theta*180./PETSC_PI),(double)angle,(double)alpha);CHKERRQ(ierr);
47786d74e61SPeter Brune   } else {
478c69d1a72SBarry Smith     ierr     = PetscInfo2(snes,"Angle %14.12e degrees exceeds threshold %14.12e, no correction applied\n",(double)(theta*180./PETSC_PI),(double)angle);CHKERRQ(ierr);
47986d74e61SPeter Brune     ierr     = VecCopy(Y,Ylast);CHKERRQ(ierr);
48086d74e61SPeter Brune     *changed = PETSC_FALSE;
481bf7f4e0aSPeter Brune   }
482bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
483bf7f4e0aSPeter Brune }
484bf7f4e0aSPeter Brune 
485bf7f4e0aSPeter Brune #undef __FUNCT__
486f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchApply"
487f40b411bSPeter Brune /*@
488cd7522eaSPeter Brune    SNESLineSearchApply - Computes the line-search update.
489f40b411bSPeter Brune 
490f1c6b773SPeter Brune    Collective on SNESLineSearch
491f40b411bSPeter Brune 
492f40b411bSPeter Brune    Input Parameters:
4938e557f58SPeter Brune +  linesearch - The linesearch context
4948e557f58SPeter Brune .  X - The current solution
4958e557f58SPeter Brune .  F - The current function
4968e557f58SPeter Brune .  fnorm - The current norm
4978e557f58SPeter Brune -  Y - The search direction
498f40b411bSPeter Brune 
499f40b411bSPeter Brune    Output Parameters:
5008e557f58SPeter Brune +  X - The new solution
5018e557f58SPeter Brune .  F - The new function
5028e557f58SPeter Brune -  fnorm - The new function norm
503f40b411bSPeter Brune 
504cd7522eaSPeter Brune    Options Database Keys:
5053c7d6663SPeter Brune + -snes_linesearch_type - basic, bt, l2, cp, shell
506cd7522eaSPeter Brune . -snes_linesearch_monitor - Print progress of line searches
5073c7d6663SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter
508cd7522eaSPeter Brune . -snes_linesearch_norms   - Turn on/off the linesearch norms
5093c7d6663SPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess
5103c7d6663SPeter Brune - -snes_linesearch_max_it - The number of iterations for iterative line searches
511cd7522eaSPeter Brune 
512cd7522eaSPeter Brune    Notes:
513cd7522eaSPeter Brune    This is typically called from within a SNESSolve() implementation in order to
514cd7522eaSPeter Brune    help with convergence of the nonlinear method.  Various SNES types use line searches
515cd7522eaSPeter Brune    in different ways, but the overarching theme is that a line search is used to determine
516cd7522eaSPeter Brune    an optimal damping parameter of a step at each iteration of the method.  Each
517cd7522eaSPeter Brune    application of the line search may invoke SNESComputeFunction several times, and
518cd7522eaSPeter Brune    therefore may be fairly expensive.
519cd7522eaSPeter Brune 
520f40b411bSPeter Brune    Level: Intermediate
521f40b411bSPeter Brune 
522f1c6b773SPeter Brune .keywords: SNESLineSearch, Create
523f40b411bSPeter Brune 
524cd7522eaSPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchPreCheck(), SNESLineSearchPostCheck(), SNESSolve(), SNESComputeFunction()
525f40b411bSPeter Brune @*/
526bf388a1fSBarry Smith PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y)
527bf388a1fSBarry Smith {
528bf7f4e0aSPeter Brune   PetscErrorCode ierr;
529bf7f4e0aSPeter Brune 
530bf388a1fSBarry Smith   PetscFunctionBegin;
531f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
532bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
533bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
534bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
535bf7f4e0aSPeter Brune 
536bf7f4e0aSPeter Brune   linesearch->success = PETSC_TRUE;
537bf7f4e0aSPeter Brune 
538bf7f4e0aSPeter Brune   linesearch->vec_sol    = X;
539bf7f4e0aSPeter Brune   linesearch->vec_update = Y;
540bf7f4e0aSPeter Brune   linesearch->vec_func   = F;
541bf7f4e0aSPeter Brune 
542f1c6b773SPeter Brune   ierr = SNESLineSearchSetUp(linesearch);CHKERRQ(ierr);
543bf7f4e0aSPeter Brune 
544f5af7f23SKarl Rupp   if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */
545bf7f4e0aSPeter Brune 
546f5af7f23SKarl Rupp   if (fnorm) linesearch->fnorm = *fnorm;
547f5af7f23SKarl Rupp   else {
548bf7f4e0aSPeter Brune     ierr = VecNorm(F, NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
549bf7f4e0aSPeter Brune   }
550bf7f4e0aSPeter Brune 
551f1c6b773SPeter Brune   ierr = PetscLogEventBegin(SNESLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
552bf7f4e0aSPeter Brune 
553bf7f4e0aSPeter Brune   ierr = (*linesearch->ops->apply)(linesearch);CHKERRQ(ierr);
554bf7f4e0aSPeter Brune 
555f1c6b773SPeter Brune   ierr = PetscLogEventEnd(SNESLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
556bf7f4e0aSPeter Brune 
557f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
558bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
559bf7f4e0aSPeter Brune }
560bf7f4e0aSPeter Brune 
561bf7f4e0aSPeter Brune #undef __FUNCT__
562f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchDestroy"
563f40b411bSPeter Brune /*@
564f1c6b773SPeter Brune    SNESLineSearchDestroy - Destroys the line search instance.
565f40b411bSPeter Brune 
566f1c6b773SPeter Brune    Collective on SNESLineSearch
567f40b411bSPeter Brune 
568f40b411bSPeter Brune    Input Parameters:
5698e557f58SPeter Brune .  linesearch - The linesearch context
570f40b411bSPeter Brune 
571f40b411bSPeter Brune    Level: Intermediate
572f40b411bSPeter Brune 
57378bcb3b5SPeter Brune .keywords: SNESLineSearch, Destroy
574f40b411bSPeter Brune 
575cd7522eaSPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchReset(), SNESDestroy()
576f40b411bSPeter Brune @*/
577bf388a1fSBarry Smith PetscErrorCode SNESLineSearchDestroy(SNESLineSearch * linesearch)
578bf388a1fSBarry Smith {
579bf7f4e0aSPeter Brune   PetscErrorCode ierr;
580bf388a1fSBarry Smith 
581bf7f4e0aSPeter Brune   PetscFunctionBegin;
582bf7f4e0aSPeter Brune   if (!*linesearch) PetscFunctionReturn(0);
583f1c6b773SPeter Brune   PetscValidHeaderSpecific((*linesearch),SNESLINESEARCH_CLASSID,1);
584bf7f4e0aSPeter Brune   if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = 0; PetscFunctionReturn(0);}
585e04113cfSBarry Smith   ierr = PetscObjectSAWsViewOff((PetscObject)*linesearch);CHKERRQ(ierr);
58622d28d08SBarry Smith   ierr = SNESLineSearchReset(*linesearch);CHKERRQ(ierr);
587f5af7f23SKarl Rupp   if ((*linesearch)->ops->destroy) (*linesearch)->ops->destroy(*linesearch);
588bf7f4e0aSPeter Brune   ierr = PetscViewerDestroy(&(*linesearch)->monitor);CHKERRQ(ierr);
589e7058c64SPeter Brune   ierr = PetscHeaderDestroy(linesearch);CHKERRQ(ierr);
590bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
591bf7f4e0aSPeter Brune }
592bf7f4e0aSPeter Brune 
593bf7f4e0aSPeter Brune #undef __FUNCT__
594f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetMonitor"
595f40b411bSPeter Brune /*@
596cd7522eaSPeter Brune    SNESLineSearchSetMonitor - Turns on/off printing useful information and debugging output about the line search.
597bf7f4e0aSPeter Brune 
598bf7f4e0aSPeter Brune    Input Parameters:
599bf7f4e0aSPeter Brune +  snes - nonlinear context obtained from SNESCreate()
600bf7f4e0aSPeter Brune -  flg - PETSC_TRUE to monitor the line search
601bf7f4e0aSPeter Brune 
602bf7f4e0aSPeter Brune    Logically Collective on SNES
603bf7f4e0aSPeter Brune 
604bf7f4e0aSPeter Brune    Options Database:
6058e557f58SPeter Brune .   -snes_linesearch_monitor - enables the monitor
606bf7f4e0aSPeter Brune 
607bf7f4e0aSPeter Brune    Level: intermediate
608bf7f4e0aSPeter Brune 
609bf7f4e0aSPeter Brune 
610cd7522eaSPeter Brune .seealso: SNESLineSearchGetMonitor(), PetscViewer
611bf7f4e0aSPeter Brune @*/
612f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetMonitor(SNESLineSearch linesearch, PetscBool flg)
613bf7f4e0aSPeter Brune {
614bf7f4e0aSPeter Brune   PetscErrorCode ierr;
615bf388a1fSBarry Smith 
616bf7f4e0aSPeter Brune   PetscFunctionBegin;
617bf7f4e0aSPeter Brune   if (flg && !linesearch->monitor) {
618ce94432eSBarry Smith     ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)linesearch),"stdout",&linesearch->monitor);CHKERRQ(ierr);
619bf7f4e0aSPeter Brune   } else if (!flg && linesearch->monitor) {
620bf7f4e0aSPeter Brune     ierr = PetscViewerDestroy(&linesearch->monitor);CHKERRQ(ierr);
621bf7f4e0aSPeter Brune   }
622bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
623bf7f4e0aSPeter Brune }
624bf7f4e0aSPeter Brune 
625bf7f4e0aSPeter Brune #undef __FUNCT__
626f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetMonitor"
627f40b411bSPeter Brune /*@
628cd7522eaSPeter Brune    SNESLineSearchGetMonitor - Gets the PetscViewer instance for the line search monitor.
6296a388c36SPeter Brune 
630f40b411bSPeter Brune    Input Parameters:
6318e557f58SPeter Brune .  linesearch - linesearch context
632f40b411bSPeter Brune 
633f40b411bSPeter Brune    Input Parameters:
6348e557f58SPeter Brune .  monitor - monitor context
635f40b411bSPeter Brune 
636f40b411bSPeter Brune    Logically Collective on SNES
637f40b411bSPeter Brune 
638f40b411bSPeter Brune 
639f40b411bSPeter Brune    Options Database Keys:
6408e557f58SPeter Brune .   -snes_linesearch_monitor - enables the monitor
641f40b411bSPeter Brune 
642f40b411bSPeter Brune    Level: intermediate
643f40b411bSPeter Brune 
644f40b411bSPeter Brune 
645cd7522eaSPeter Brune .seealso: SNESLineSearchSetMonitor(), PetscViewer
646f40b411bSPeter Brune @*/
647f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetMonitor(SNESLineSearch linesearch, PetscViewer *monitor)
6486a388c36SPeter Brune {
6496a388c36SPeter Brune   PetscFunctionBegin;
650f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
6516a388c36SPeter Brune   if (monitor) {
6526a388c36SPeter Brune     PetscValidPointer(monitor, 2);
6536a388c36SPeter Brune     *monitor = linesearch->monitor;
6546a388c36SPeter Brune   }
6556a388c36SPeter Brune   PetscFunctionReturn(0);
6566a388c36SPeter Brune }
6576a388c36SPeter Brune 
6586a388c36SPeter Brune #undef __FUNCT__
659f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetFromOptions"
660f40b411bSPeter Brune /*@
661f1c6b773SPeter Brune    SNESLineSearchSetFromOptions - Sets options for the line search
662f40b411bSPeter Brune 
663f40b411bSPeter Brune    Input Parameters:
6648e557f58SPeter Brune .  linesearch - linesearch context
665f40b411bSPeter Brune 
666f40b411bSPeter Brune    Options Database Keys:
6673c7d6663SPeter Brune + -snes_linesearch_type <type> - basic, bt, l2, cp, shell
6683c7d6663SPeter Brune . -snes_linesearch_order <order> - 1, 2, 3.  Most types only support certain orders (bt supports 2 or 3)
6695a9b6599SPeter Brune . -snes_linesearch_norms   - Turn on/off the linesearch norms for the basic linesearch type
67071eef1aeSPeter Brune . -snes_linesearch_minlambda - The minimum step length
6711a9b3a06SPeter Brune . -snes_linesearch_maxstep - The maximum step size
6721a9b3a06SPeter Brune . -snes_linesearch_rtol - Relative tolerance for iterative line searches
6731a9b3a06SPeter Brune . -snes_linesearch_atol - Absolute tolerance for iterative line searches
6741a9b3a06SPeter Brune . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches
6751a9b3a06SPeter Brune . -snes_linesearch_max_it - The number of iterations for iterative line searches
676cd7522eaSPeter Brune . -snes_linesearch_monitor - Print progress of line searches
6778e557f58SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter
678cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess.
6791a9b3a06SPeter Brune . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method
6801a9b3a06SPeter Brune - -snes_linesearch_precheck_picard_angle - Angle used in picard precheck method
681f40b411bSPeter Brune 
682f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
683f40b411bSPeter Brune 
684f40b411bSPeter Brune    Level: intermediate
685f40b411bSPeter Brune 
6863c7d6663SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetOrder(), SNESLineSearchSetType(), SNESLineSearchSetTolerances(), SNESLineSearchSetDamping(), SNESLineSearchPreCheckPicard()
687f40b411bSPeter Brune @*/
688bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch)
689bf388a1fSBarry Smith {
690bf7f4e0aSPeter Brune   PetscErrorCode ierr;
6911a4f838cSPeter Brune   const char     *deft = SNESLINESEARCHBASIC;
692bf7f4e0aSPeter Brune   char           type[256];
693bf7f4e0aSPeter Brune   PetscBool      flg, set;
694bf388a1fSBarry Smith 
695bf7f4e0aSPeter Brune   PetscFunctionBegin;
696607a6623SBarry Smith   if (!SNESLineSearchRegisterAllCalled) {ierr = SNESLineSearchRegisterAll();CHKERRQ(ierr);}
697bf7f4e0aSPeter Brune 
698bf7f4e0aSPeter Brune   ierr = PetscObjectOptionsBegin((PetscObject)linesearch);CHKERRQ(ierr);
699f5af7f23SKarl Rupp   if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name;
700a264d7a6SBarry Smith   ierr = PetscOptionsFList("-snes_linesearch_type","Linesearch type","SNESLineSearchSetType",SNESLineSearchList,deft,type,256,&flg);CHKERRQ(ierr);
701bf7f4e0aSPeter Brune   if (flg) {
702f1c6b773SPeter Brune     ierr = SNESLineSearchSetType(linesearch,type);CHKERRQ(ierr);
703bf7f4e0aSPeter Brune   } else if (!((PetscObject)linesearch)->type_name) {
704f1c6b773SPeter Brune     ierr = SNESLineSearchSetType(linesearch,deft);CHKERRQ(ierr);
705bf7f4e0aSPeter Brune   }
706bf7f4e0aSPeter Brune 
7077a35526eSPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_monitor","Print progress of line searches","SNESSNESLineSearchSetMonitor",
708bf7f4e0aSPeter Brune                           linesearch->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
709f1c6b773SPeter Brune   if (set) {ierr = SNESLineSearchSetMonitor(linesearch,flg);CHKERRQ(ierr);}
710bf7f4e0aSPeter Brune 
7111a9b3a06SPeter Brune   /* tolerances */
71271eef1aeSPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_minlambda","Minimum step length","SNESLineSearchSetTolerances",linesearch->steptol,&linesearch->steptol,0);CHKERRQ(ierr);
7131a9b3a06SPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_maxstep","Maximum step size","SNESLineSearchSetTolerances",linesearch->maxstep,&linesearch->maxstep,0);CHKERRQ(ierr);
7141a9b3a06SPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_rtol","Relative tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->rtol,&linesearch->rtol,0);CHKERRQ(ierr);
7151a9b3a06SPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_atol","Absolute tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->atol,&linesearch->atol,0);CHKERRQ(ierr);
7161a9b3a06SPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_ltol","Change in lambda tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->ltol,&linesearch->ltol,0);CHKERRQ(ierr);
7178e557f58SPeter Brune   ierr = PetscOptionsInt("-snes_linesearch_max_it","Maximum iterations for iterative line searches","SNESLineSearchSetTolerances",linesearch->max_its,&linesearch->max_its,0);CHKERRQ(ierr);
7187a35526eSPeter Brune 
7191a9b3a06SPeter Brune   /* damping parameters */
7201a9b3a06SPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_damping","Line search damping and initial step guess","SNESLineSearchSetDamping",linesearch->damping,&linesearch->damping,0);CHKERRQ(ierr);
7211a9b3a06SPeter Brune 
7221a9b3a06SPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_keeplambda","Use previous lambda as damping","SNESLineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,0);CHKERRQ(ierr);
7231a9b3a06SPeter Brune 
7241a9b3a06SPeter Brune   /* precheck */
7257a35526eSPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_precheck_picard","Use a correction that sometimes improves convergence of Picard iteration","SNESLineSearchPreCheckPicard",flg,&flg,&set);CHKERRQ(ierr);
7267a35526eSPeter Brune   if (set) {
7277a35526eSPeter Brune     if (flg) {
7287a35526eSPeter Brune       linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */
729f5af7f23SKarl Rupp 
7307a35526eSPeter Brune       ierr = PetscOptionsReal("-snes_linesearch_precheck_picard_angle","Maximum angle at which to activate the correction",
7310298fd71SBarry Smith                               "none",linesearch->precheck_picard_angle,&linesearch->precheck_picard_angle,NULL);CHKERRQ(ierr);
7327a35526eSPeter Brune       ierr = SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&linesearch->precheck_picard_angle);CHKERRQ(ierr);
7337a35526eSPeter Brune     } else {
7340298fd71SBarry Smith       ierr = SNESLineSearchSetPreCheck(linesearch,NULL,NULL);CHKERRQ(ierr);
7357a35526eSPeter Brune     }
7367a35526eSPeter Brune   }
737b000cd8dSPeter Brune   ierr = PetscOptionsInt("-snes_linesearch_order","Order of approximation used in the line search","SNESLineSearchSetOrder",linesearch->order,&linesearch->order,0);CHKERRQ(ierr);
7381a9b3a06SPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_norms","Compute final norms in line search","SNESLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,0);CHKERRQ(ierr);
7397a35526eSPeter Brune 
7405a9b6599SPeter Brune   if (linesearch->ops->setfromoptions) {
7415a9b6599SPeter Brune     (*linesearch->ops->setfromoptions)(linesearch);CHKERRQ(ierr);
7425a9b6599SPeter Brune   }
7435a9b6599SPeter Brune 
744bf7f4e0aSPeter Brune   ierr = PetscObjectProcessOptionsHandlers((PetscObject)linesearch);CHKERRQ(ierr);
745bf7f4e0aSPeter Brune   ierr = PetscOptionsEnd();CHKERRQ(ierr);
746bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
747bf7f4e0aSPeter Brune }
748bf7f4e0aSPeter Brune 
749bf7f4e0aSPeter Brune #undef __FUNCT__
750f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchView"
751f40b411bSPeter Brune /*@
752cd7522eaSPeter Brune    SNESLineSearchView - Prints useful information about the line search not
753cd7522eaSPeter Brune    related to an individual call.
754f40b411bSPeter Brune 
755f40b411bSPeter Brune    Input Parameters:
7568e557f58SPeter Brune .  linesearch - linesearch context
757f40b411bSPeter Brune 
758f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
759f40b411bSPeter Brune 
760f40b411bSPeter Brune    Level: intermediate
761f40b411bSPeter Brune 
762f1c6b773SPeter Brune .seealso: SNESLineSearchCreate()
763f40b411bSPeter Brune @*/
764bf388a1fSBarry Smith PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer)
765bf388a1fSBarry Smith {
7667f1410a3SPeter Brune   PetscErrorCode ierr;
7677f1410a3SPeter Brune   PetscBool      iascii;
768bf388a1fSBarry Smith 
769bf7f4e0aSPeter Brune   PetscFunctionBegin;
7707f1410a3SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
7717f1410a3SPeter Brune   if (!viewer) {
772ce94432eSBarry Smith     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch),&viewer);CHKERRQ(ierr);
7737f1410a3SPeter Brune   }
7747f1410a3SPeter Brune   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
7757f1410a3SPeter Brune   PetscCheckSameComm(linesearch,1,viewer,2);
776f40b411bSPeter Brune 
777251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
7787f1410a3SPeter Brune   if (iascii) {
779dae58748SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)linesearch,viewer);CHKERRQ(ierr);
7807f1410a3SPeter Brune     if (linesearch->ops->view) {
7817f1410a3SPeter Brune       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
7827f1410a3SPeter Brune       ierr = (*linesearch->ops->view)(linesearch,viewer);CHKERRQ(ierr);
7837f1410a3SPeter Brune       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
7847f1410a3SPeter Brune     }
785c69d1a72SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep,(double)linesearch->steptol);CHKERRQ(ierr);
786c69d1a72SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol,(double)linesearch->atol,(double)linesearch->ltol);CHKERRQ(ierr);
7877f1410a3SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D\n", linesearch->max_its);CHKERRQ(ierr);
7886b2b7091SBarry Smith     if (linesearch->ops->precheck) {
7896b2b7091SBarry Smith       if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) {
7907f1410a3SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  using precheck step to speed up Picard convergence\n", linesearch->max_its);CHKERRQ(ierr);
7917f1410a3SPeter Brune       } else {
7927f1410a3SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  using user-defined precheck step\n", linesearch->max_its);CHKERRQ(ierr);
7937f1410a3SPeter Brune       }
7947f1410a3SPeter Brune     }
7956b2b7091SBarry Smith     if (linesearch->ops->postcheck) {
7967f1410a3SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer,"  using user-defined postcheck step\n", linesearch->max_its);CHKERRQ(ierr);
7977f1410a3SPeter Brune     }
7987f1410a3SPeter Brune   }
799bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
800bf7f4e0aSPeter Brune }
801bf7f4e0aSPeter Brune 
802bf7f4e0aSPeter Brune #undef __FUNCT__
803f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetType"
804ea5d4fccSPeter Brune /*@C
805f1c6b773SPeter Brune    SNESLineSearchSetType - Sets the linesearch type
806f40b411bSPeter Brune 
807f40b411bSPeter Brune    Input Parameters:
8088e557f58SPeter Brune +  linesearch - linesearch context
809f40b411bSPeter Brune -  type - The type of line search to be used
810f40b411bSPeter Brune 
811cd7522eaSPeter Brune    Available Types:
812cd7522eaSPeter Brune +  basic - Simple damping line search.
8138e557f58SPeter Brune .  bt - Backtracking line search over the L2 norm of the function
8148e557f58SPeter Brune .  l2 - Secant line search over the L2 norm of the function
8158e557f58SPeter Brune .  cp - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x)
8168e557f58SPeter Brune -  shell - User provided SNESLineSearch implementation
817cd7522eaSPeter Brune 
818f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
819f40b411bSPeter Brune 
820f40b411bSPeter Brune    Level: intermediate
821f40b411bSPeter Brune 
822f40b411bSPeter Brune 
823f1c6b773SPeter Brune .seealso: SNESLineSearchCreate()
824f40b411bSPeter Brune @*/
82519fd82e9SBarry Smith PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type)
826bf7f4e0aSPeter Brune {
827f1c6b773SPeter Brune   PetscErrorCode ierr,(*r)(SNESLineSearch);
828bf7f4e0aSPeter Brune   PetscBool      match;
829bf7f4e0aSPeter Brune 
830bf7f4e0aSPeter Brune   PetscFunctionBegin;
831f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
832bf7f4e0aSPeter Brune   PetscValidCharPointer(type,2);
833bf7f4e0aSPeter Brune 
834251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)linesearch,type,&match);CHKERRQ(ierr);
835bf7f4e0aSPeter Brune   if (match) PetscFunctionReturn(0);
836bf7f4e0aSPeter Brune 
8371c9cd337SJed Brown   ierr = PetscFunctionListFind(SNESLineSearchList,type,&r);CHKERRQ(ierr);
838bf7f4e0aSPeter Brune   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type);
839bf7f4e0aSPeter Brune   /* Destroy the previous private linesearch context */
840bf7f4e0aSPeter Brune   if (linesearch->ops->destroy) {
841bf7f4e0aSPeter Brune     ierr = (*(linesearch)->ops->destroy)(linesearch);CHKERRQ(ierr);
842f5af7f23SKarl Rupp 
8430298fd71SBarry Smith     linesearch->ops->destroy = NULL;
844bf7f4e0aSPeter Brune   }
845f1c6b773SPeter Brune   /* Reinitialize function pointers in SNESLineSearchOps structure */
846bf7f4e0aSPeter Brune   linesearch->ops->apply          = 0;
847bf7f4e0aSPeter Brune   linesearch->ops->view           = 0;
848bf7f4e0aSPeter Brune   linesearch->ops->setfromoptions = 0;
849bf7f4e0aSPeter Brune   linesearch->ops->destroy        = 0;
850bf7f4e0aSPeter Brune 
851bf7f4e0aSPeter Brune   ierr = PetscObjectChangeTypeName((PetscObject)linesearch,type);CHKERRQ(ierr);
852bf7f4e0aSPeter Brune   ierr = (*r)(linesearch);CHKERRQ(ierr);
853bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
854bf7f4e0aSPeter Brune }
855bf7f4e0aSPeter Brune 
856bf7f4e0aSPeter Brune #undef __FUNCT__
857f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetSNES"
858f40b411bSPeter Brune /*@
85978bcb3b5SPeter Brune    SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation.
860f40b411bSPeter Brune 
861f40b411bSPeter Brune    Input Parameters:
8628e557f58SPeter Brune +  linesearch - linesearch context
863f40b411bSPeter Brune -  snes - The snes instance
864f40b411bSPeter Brune 
86578bcb3b5SPeter Brune    Level: developer
86678bcb3b5SPeter Brune 
86778bcb3b5SPeter Brune    Notes:
86878bcb3b5SPeter Brune    This happens automatically when the line search is gotten/created with
8697601faf0SJed Brown    SNESGetLineSearch().  This routine is therefore mainly called within SNES
87078bcb3b5SPeter Brune    implementations.
871f40b411bSPeter Brune 
8728141a3b9SPeter Brune    Level: developer
873f40b411bSPeter Brune 
874cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
875f40b411bSPeter Brune @*/
876bf388a1fSBarry Smith PetscErrorCode  SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes)
877bf388a1fSBarry Smith {
878bf7f4e0aSPeter Brune   PetscFunctionBegin;
879f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
880bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,2);
881bf7f4e0aSPeter Brune   linesearch->snes = snes;
882bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
883bf7f4e0aSPeter Brune }
884bf7f4e0aSPeter Brune 
885bf7f4e0aSPeter Brune #undef __FUNCT__
886f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetSNES"
887f40b411bSPeter Brune /*@
8888141a3b9SPeter Brune    SNESLineSearchGetSNES - Gets the SNES instance associated with the line search.
8898141a3b9SPeter Brune    Having an associated SNES is necessary because most line search implementations must be able to
8908141a3b9SPeter Brune    evaluate the function using SNESComputeFunction() for the associated SNES.  This routine
8918141a3b9SPeter Brune    is used in the line search implementations when one must get this associated SNES instance.
892f40b411bSPeter Brune 
893f40b411bSPeter Brune    Input Parameters:
8948e557f58SPeter Brune .  linesearch - linesearch context
895f40b411bSPeter Brune 
896f40b411bSPeter Brune    Output Parameters:
897f40b411bSPeter Brune .  snes - The snes instance
898f40b411bSPeter Brune 
8998141a3b9SPeter Brune    Level: developer
900f40b411bSPeter Brune 
901cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
902f40b411bSPeter Brune @*/
903bf388a1fSBarry Smith PetscErrorCode  SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes)
904bf388a1fSBarry Smith {
905bf7f4e0aSPeter Brune   PetscFunctionBegin;
906f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
9076a388c36SPeter Brune   PetscValidPointer(snes, 2);
908bf7f4e0aSPeter Brune   *snes = linesearch->snes;
909bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
910bf7f4e0aSPeter Brune }
911bf7f4e0aSPeter Brune 
9126a388c36SPeter Brune #undef __FUNCT__
913f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetLambda"
914f40b411bSPeter Brune /*@
915f1c6b773SPeter Brune    SNESLineSearchGetLambda - Gets the last linesearch steplength discovered.
916f40b411bSPeter Brune 
917f40b411bSPeter Brune    Input Parameters:
9188e557f58SPeter Brune .  linesearch - linesearch context
919f40b411bSPeter Brune 
920f40b411bSPeter Brune    Output Parameters:
921cd7522eaSPeter Brune .  lambda - The last steplength computed during SNESLineSearchApply()
922f40b411bSPeter Brune 
92378bcb3b5SPeter Brune    Level: advanced
92478bcb3b5SPeter Brune 
9258e557f58SPeter Brune    Notes:
9268e557f58SPeter Brune    This is useful in methods where the solver is ill-scaled and
92778bcb3b5SPeter Brune    requires some adaptive notion of the difference in scale between the
92878bcb3b5SPeter Brune    solution and the function.  For instance, SNESQN may be scaled by the
92978bcb3b5SPeter Brune    line search lambda using the argument -snes_qn_scaling ls.
93078bcb3b5SPeter Brune 
931f40b411bSPeter Brune 
932cd7522eaSPeter Brune .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetDamping(), SNESLineSearchApply()
933f40b411bSPeter Brune @*/
934f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda)
9356a388c36SPeter Brune {
9366a388c36SPeter Brune   PetscFunctionBegin;
937f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
9386a388c36SPeter Brune   PetscValidPointer(lambda, 2);
9396a388c36SPeter Brune   *lambda = linesearch->lambda;
9406a388c36SPeter Brune   PetscFunctionReturn(0);
9416a388c36SPeter Brune }
9426a388c36SPeter Brune 
9436a388c36SPeter Brune #undef __FUNCT__
944f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetLambda"
945f40b411bSPeter Brune /*@
946f1c6b773SPeter Brune    SNESLineSearchSetLambda - Sets the linesearch steplength.
947f40b411bSPeter Brune 
948f40b411bSPeter Brune    Input Parameters:
9498e557f58SPeter Brune +  linesearch - linesearch context
950f40b411bSPeter Brune -  lambda - The last steplength.
951f40b411bSPeter Brune 
952cd7522eaSPeter Brune    Notes:
953cd7522eaSPeter Brune    This routine is typically used within implementations of SNESLineSearchApply
954cd7522eaSPeter Brune    to set the final steplength.  This routine (and SNESLineSearchGetLambda()) were
955cd7522eaSPeter Brune    added in order to facilitate Quasi-Newton methods that use the previous steplength
956cd7522eaSPeter Brune    as an inner scaling parameter.
957cd7522eaSPeter Brune 
95878bcb3b5SPeter Brune    Level: advanced
959f40b411bSPeter Brune 
960f1c6b773SPeter Brune .seealso: SNESLineSearchGetLambda()
961f40b411bSPeter Brune @*/
962f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
9636a388c36SPeter Brune {
9646a388c36SPeter Brune   PetscFunctionBegin;
965f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
9666a388c36SPeter Brune   linesearch->lambda = lambda;
9676a388c36SPeter Brune   PetscFunctionReturn(0);
9686a388c36SPeter Brune }
9696a388c36SPeter Brune 
9706a388c36SPeter Brune #undef  __FUNCT__
971f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetTolerances"
972f40b411bSPeter Brune /*@
9733c7d6663SPeter Brune    SNESLineSearchGetTolerances - Gets the tolerances for the linesearch.  These include
97478bcb3b5SPeter Brune    tolerances for the relative and absolute change in the function norm, the change
97578bcb3b5SPeter Brune    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
97678bcb3b5SPeter Brune    and the maximum number of iterations the line search procedure may take.
977f40b411bSPeter Brune 
978f40b411bSPeter Brune    Input Parameters:
9798e557f58SPeter Brune .  linesearch - linesearch context
980f40b411bSPeter Brune 
981f40b411bSPeter Brune    Output Parameters:
982516fe3c3SPeter Brune +  steptol - The minimum steplength
9836cc8e53bSPeter Brune .  maxstep - The maximum steplength
984516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
985516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
986516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
987516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
988f40b411bSPeter Brune 
98978bcb3b5SPeter Brune    Level: intermediate
99078bcb3b5SPeter Brune 
99178bcb3b5SPeter Brune    Notes:
99278bcb3b5SPeter Brune    Different line searches may implement these parameters slightly differently as
9933c7d6663SPeter Brune    the type requires.
994516fe3c3SPeter Brune 
995f1c6b773SPeter Brune .seealso: SNESLineSearchSetTolerances()
996f40b411bSPeter Brune @*/
997f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
9986a388c36SPeter Brune {
9996a388c36SPeter Brune   PetscFunctionBegin;
1000f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1001516fe3c3SPeter Brune   if (steptol) {
10026a388c36SPeter Brune     PetscValidPointer(steptol, 2);
10036a388c36SPeter Brune     *steptol = linesearch->steptol;
1004516fe3c3SPeter Brune   }
1005516fe3c3SPeter Brune   if (maxstep) {
1006516fe3c3SPeter Brune     PetscValidPointer(maxstep, 3);
1007516fe3c3SPeter Brune     *maxstep = linesearch->maxstep;
1008516fe3c3SPeter Brune   }
1009516fe3c3SPeter Brune   if (rtol) {
1010516fe3c3SPeter Brune     PetscValidPointer(rtol, 4);
1011516fe3c3SPeter Brune     *rtol = linesearch->rtol;
1012516fe3c3SPeter Brune   }
1013516fe3c3SPeter Brune   if (atol) {
1014516fe3c3SPeter Brune     PetscValidPointer(atol, 5);
1015516fe3c3SPeter Brune     *atol = linesearch->atol;
1016516fe3c3SPeter Brune   }
1017516fe3c3SPeter Brune   if (ltol) {
1018516fe3c3SPeter Brune     PetscValidPointer(ltol, 6);
1019516fe3c3SPeter Brune     *ltol = linesearch->ltol;
1020516fe3c3SPeter Brune   }
1021516fe3c3SPeter Brune   if (max_its) {
1022516fe3c3SPeter Brune     PetscValidPointer(max_its, 7);
1023516fe3c3SPeter Brune     *max_its = linesearch->max_its;
1024516fe3c3SPeter Brune   }
10256a388c36SPeter Brune   PetscFunctionReturn(0);
10266a388c36SPeter Brune }
10276a388c36SPeter Brune 
10286a388c36SPeter Brune #undef  __FUNCT__
1029f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetTolerances"
1030f40b411bSPeter Brune /*@
10313c7d6663SPeter Brune    SNESLineSearchSetTolerances -  Gets the tolerances for the linesearch.  These include
103278bcb3b5SPeter Brune    tolerances for the relative and absolute change in the function norm, the change
103378bcb3b5SPeter Brune    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
103478bcb3b5SPeter Brune    and the maximum number of iterations the line search procedure may take.
1035f40b411bSPeter Brune 
1036f40b411bSPeter Brune    Input Parameters:
10378e557f58SPeter Brune +  linesearch - linesearch context
1038516fe3c3SPeter Brune .  steptol - The minimum steplength
10396cc8e53bSPeter Brune .  maxstep - The maximum steplength
1040516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
1041516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
1042516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
1043516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
1044f40b411bSPeter Brune 
104578bcb3b5SPeter Brune    Notes:
10463c7d6663SPeter Brune    The user may choose to not set any of the tolerances using PETSC_DEFAULT in
104778bcb3b5SPeter Brune    place of an argument.
1048f40b411bSPeter Brune 
104978bcb3b5SPeter Brune    Level: intermediate
1050516fe3c3SPeter Brune 
1051f1c6b773SPeter Brune .seealso: SNESLineSearchGetTolerances()
1052f40b411bSPeter Brune @*/
1053f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its)
10546a388c36SPeter Brune {
10556a388c36SPeter Brune   PetscFunctionBegin;
1056f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1057d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,steptol,2);
1058d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,maxstep,3);
1059d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,rtol,4);
1060d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,atol,5);
1061d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,ltol,6);
1062d3952184SSatish Balay   PetscValidLogicalCollectiveInt(linesearch,max_its,7);
1063d3952184SSatish Balay 
1064d3952184SSatish Balay   if (steptol!= PETSC_DEFAULT) {
1065ce94432eSBarry Smith     if (steptol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Minimum step length %14.12e must be non-negative",(double)steptol);
10666a388c36SPeter Brune     linesearch->steptol = steptol;
1067d3952184SSatish Balay   }
1068d3952184SSatish Balay 
1069d3952184SSatish Balay   if (maxstep!= PETSC_DEFAULT) {
1070ce94432eSBarry Smith     if (maxstep < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum step length %14.12e must be non-negative",(double)maxstep);
1071516fe3c3SPeter Brune     linesearch->maxstep = maxstep;
1072d3952184SSatish Balay   }
1073d3952184SSatish Balay 
1074d3952184SSatish Balay   if (rtol != PETSC_DEFAULT) {
1075ce94432eSBarry Smith     if (rtol < 0.0 || 1.0 <= rtol) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Relative tolerance %14.12e must be non-negative and less than 1.0",(double)rtol);
1076516fe3c3SPeter Brune     linesearch->rtol = rtol;
1077d3952184SSatish Balay   }
1078d3952184SSatish Balay 
1079d3952184SSatish Balay   if (atol != PETSC_DEFAULT) {
1080ce94432eSBarry Smith     if (atol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %14.12e must be non-negative",(double)atol);
1081516fe3c3SPeter Brune     linesearch->atol = atol;
1082d3952184SSatish Balay   }
1083d3952184SSatish Balay 
1084d3952184SSatish Balay   if (ltol != PETSC_DEFAULT) {
1085ce94432eSBarry Smith     if (ltol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Labmda tolerance %14.12e must be non-negative",(double)ltol);
1086516fe3c3SPeter Brune     linesearch->ltol = ltol;
1087d3952184SSatish Balay   }
1088d3952184SSatish Balay 
1089d3952184SSatish Balay   if (max_its != PETSC_DEFAULT) {
1090ce94432eSBarry Smith     if (max_its < 0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",max_its);
1091516fe3c3SPeter Brune     linesearch->max_its = max_its;
1092d3952184SSatish Balay   }
10936a388c36SPeter Brune   PetscFunctionReturn(0);
10946a388c36SPeter Brune }
10956a388c36SPeter Brune 
10966a388c36SPeter Brune #undef __FUNCT__
1097f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetDamping"
1098f40b411bSPeter Brune /*@
1099f1c6b773SPeter Brune    SNESLineSearchGetDamping - Gets the line search damping parameter.
1100f40b411bSPeter Brune 
1101f40b411bSPeter Brune    Input Parameters:
11028e557f58SPeter Brune .  linesearch - linesearch context
1103f40b411bSPeter Brune 
1104f40b411bSPeter Brune    Output Parameters:
11058e557f58SPeter Brune .  damping - The damping parameter
1106f40b411bSPeter Brune 
110778bcb3b5SPeter Brune    Level: advanced
1108f40b411bSPeter Brune 
110978bcb3b5SPeter Brune .seealso: SNESLineSearchGetStepTolerance(), SNESQN
1110f40b411bSPeter Brune @*/
1111f40b411bSPeter Brune 
1112f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping)
11136a388c36SPeter Brune {
11146a388c36SPeter Brune   PetscFunctionBegin;
1115f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
11166a388c36SPeter Brune   PetscValidPointer(damping, 2);
11176a388c36SPeter Brune   *damping = linesearch->damping;
11186a388c36SPeter Brune   PetscFunctionReturn(0);
11196a388c36SPeter Brune }
11206a388c36SPeter Brune 
11216a388c36SPeter Brune #undef __FUNCT__
1122f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetDamping"
1123f40b411bSPeter Brune /*@
1124f1c6b773SPeter Brune    SNESLineSearchSetDamping - Sets the line search damping paramter.
1125f40b411bSPeter Brune 
1126f40b411bSPeter Brune    Input Parameters:
112778bcb3b5SPeter Brune .  linesearch - linesearch context
112878bcb3b5SPeter Brune .  damping - The damping parameter
1129f40b411bSPeter Brune 
1130f40b411bSPeter Brune    Level: intermediate
1131f40b411bSPeter Brune 
1132cd7522eaSPeter Brune    Notes:
1133cd7522eaSPeter Brune    The basic line search merely takes the update step scaled by the damping parameter.
1134cd7522eaSPeter Brune    The use of the damping parameter in the l2 and cp line searches is much more subtle;
113578bcb3b5SPeter Brune    it is used as a starting point in calculating the secant step. However, the eventual
1136cd7522eaSPeter Brune    step may be of greater length than the damping parameter.  In the bt line search it is
1137cd7522eaSPeter Brune    used as the maximum possible step length, as the bt line search only backtracks.
1138cd7522eaSPeter Brune 
1139f1c6b773SPeter Brune .seealso: SNESLineSearchGetDamping()
1140f40b411bSPeter Brune @*/
1141f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping)
11426a388c36SPeter Brune {
11436a388c36SPeter Brune   PetscFunctionBegin;
1144f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
11456a388c36SPeter Brune   linesearch->damping = damping;
11466a388c36SPeter Brune   PetscFunctionReturn(0);
11476a388c36SPeter Brune }
11486a388c36SPeter Brune 
11496a388c36SPeter Brune #undef __FUNCT__
115059405d5eSPeter Brune #define __FUNCT__ "SNESLineSearchGetOrder"
115159405d5eSPeter Brune /*@
115259405d5eSPeter Brune    SNESLineSearchGetOrder - Gets the line search approximation order.
115359405d5eSPeter Brune 
115459405d5eSPeter Brune    Input Parameters:
115578bcb3b5SPeter Brune .  linesearch - linesearch context
115659405d5eSPeter Brune 
115759405d5eSPeter Brune    Output Parameters:
11588e557f58SPeter Brune .  order - The order
115959405d5eSPeter Brune 
116078bcb3b5SPeter Brune    Possible Values for order:
11613c7d6663SPeter Brune +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
11623c7d6663SPeter Brune .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
11633c7d6663SPeter Brune -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
116478bcb3b5SPeter Brune 
116559405d5eSPeter Brune    Level: intermediate
116659405d5eSPeter Brune 
116759405d5eSPeter Brune .seealso: SNESLineSearchSetOrder()
116859405d5eSPeter Brune @*/
116959405d5eSPeter Brune 
1170b000cd8dSPeter Brune PetscErrorCode  SNESLineSearchGetOrder(SNESLineSearch linesearch,PetscInt *order)
117159405d5eSPeter Brune {
117259405d5eSPeter Brune   PetscFunctionBegin;
117359405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
117459405d5eSPeter Brune   PetscValidPointer(order, 2);
117559405d5eSPeter Brune   *order = linesearch->order;
117659405d5eSPeter Brune   PetscFunctionReturn(0);
117759405d5eSPeter Brune }
117859405d5eSPeter Brune 
117959405d5eSPeter Brune #undef __FUNCT__
118059405d5eSPeter Brune #define __FUNCT__ "SNESLineSearchSetOrder"
118159405d5eSPeter Brune /*@
118259405d5eSPeter Brune    SNESLineSearchSetOrder - Sets the line search damping paramter.
118359405d5eSPeter Brune 
118459405d5eSPeter Brune    Input Parameters:
118578bcb3b5SPeter Brune .  linesearch - linesearch context
118678bcb3b5SPeter Brune .  order - The damping parameter
118759405d5eSPeter Brune 
118859405d5eSPeter Brune    Level: intermediate
118959405d5eSPeter Brune 
119078bcb3b5SPeter Brune    Possible Values for order:
11913c7d6663SPeter Brune +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
11923c7d6663SPeter Brune .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
11933c7d6663SPeter Brune -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
119478bcb3b5SPeter Brune 
119559405d5eSPeter Brune    Notes:
119659405d5eSPeter Brune    Variable orders are supported by the following line searches:
119778bcb3b5SPeter Brune +  bt - cubic and quadratic
119878bcb3b5SPeter Brune -  cp - linear and quadratic
119959405d5eSPeter Brune 
120059405d5eSPeter Brune .seealso: SNESLineSearchGetOrder()
120159405d5eSPeter Brune @*/
1202b000cd8dSPeter Brune PetscErrorCode  SNESLineSearchSetOrder(SNESLineSearch linesearch,PetscInt order)
120359405d5eSPeter Brune {
120459405d5eSPeter Brune   PetscFunctionBegin;
120559405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
120659405d5eSPeter Brune   linesearch->order = order;
120759405d5eSPeter Brune   PetscFunctionReturn(0);
120859405d5eSPeter Brune }
120959405d5eSPeter Brune 
121059405d5eSPeter Brune #undef __FUNCT__
1211f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetNorms"
1212f40b411bSPeter Brune /*@
1213f1c6b773SPeter Brune    SNESLineSearchGetNorms - Gets the norms for for X, Y, and F.
1214f40b411bSPeter Brune 
1215f40b411bSPeter Brune    Input Parameters:
121678bcb3b5SPeter Brune .  linesearch - linesearch context
1217f40b411bSPeter Brune 
1218f40b411bSPeter Brune    Output Parameters:
1219f40b411bSPeter Brune +  xnorm - The norm of the current solution
1220f40b411bSPeter Brune .  fnorm - The norm of the current function
1221f40b411bSPeter Brune -  ynorm - The norm of the current update
1222f40b411bSPeter Brune 
1223cd7522eaSPeter Brune    Notes:
1224cd7522eaSPeter Brune    This function is mainly called from SNES implementations.
1225cd7522eaSPeter Brune 
122678bcb3b5SPeter Brune    Level: developer
1227f40b411bSPeter Brune 
1228f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms() SNESLineSearchGetVecs()
1229f40b411bSPeter Brune @*/
1230f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm)
1231bf7f4e0aSPeter Brune {
1232bf7f4e0aSPeter Brune   PetscFunctionBegin;
1233f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1234f5af7f23SKarl Rupp   if (xnorm) *xnorm = linesearch->xnorm;
1235f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
1236f5af7f23SKarl Rupp   if (ynorm) *ynorm = linesearch->ynorm;
1237bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1238bf7f4e0aSPeter Brune }
1239bf7f4e0aSPeter Brune 
1240e7058c64SPeter Brune #undef __FUNCT__
1241f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetNorms"
1242f40b411bSPeter Brune /*@
1243f1c6b773SPeter Brune    SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F.
1244f40b411bSPeter Brune 
1245f40b411bSPeter Brune    Input Parameters:
124678bcb3b5SPeter Brune +  linesearch - linesearch context
1247f40b411bSPeter Brune .  xnorm - The norm of the current solution
1248f40b411bSPeter Brune .  fnorm - The norm of the current function
1249f40b411bSPeter Brune -  ynorm - The norm of the current update
1250f40b411bSPeter Brune 
125178bcb3b5SPeter Brune    Level: advanced
1252f40b411bSPeter Brune 
1253f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1254f40b411bSPeter Brune @*/
1255f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
12566a388c36SPeter Brune {
12576a388c36SPeter Brune   PetscFunctionBegin;
1258f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
12596a388c36SPeter Brune   linesearch->xnorm = xnorm;
12606a388c36SPeter Brune   linesearch->fnorm = fnorm;
12616a388c36SPeter Brune   linesearch->ynorm = ynorm;
12626a388c36SPeter Brune   PetscFunctionReturn(0);
12636a388c36SPeter Brune }
12646a388c36SPeter Brune 
12656a388c36SPeter Brune #undef __FUNCT__
1266f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchComputeNorms"
1267f40b411bSPeter Brune /*@
1268f1c6b773SPeter Brune    SNESLineSearchComputeNorms - Computes the norms of X, F, and Y.
1269f40b411bSPeter Brune 
1270f40b411bSPeter Brune    Input Parameters:
127178bcb3b5SPeter Brune .  linesearch - linesearch context
1272f40b411bSPeter Brune 
1273f40b411bSPeter Brune    Options Database Keys:
12748e557f58SPeter Brune .   -snes_linesearch_norms - turn norm computation on or off
1275f40b411bSPeter Brune 
1276f40b411bSPeter Brune    Level: intermediate
1277f40b411bSPeter Brune 
127878bcb3b5SPeter Brune .seealso: SNESLineSearchGetNorms, SNESLineSearchSetNorms(), SNESLineSearchSetComputeNorms()
1279f40b411bSPeter Brune @*/
1280f1c6b773SPeter Brune PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
12816a388c36SPeter Brune {
12826a388c36SPeter Brune   PetscErrorCode ierr;
12839bd66eb0SPeter Brune   SNES           snes;
1284bf388a1fSBarry Smith 
12856a388c36SPeter Brune   PetscFunctionBegin;
12866a388c36SPeter Brune   if (linesearch->norms) {
12879bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
1288f1c6b773SPeter Brune       ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
12899bd66eb0SPeter Brune       ierr = VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
12909bd66eb0SPeter Brune       ierr = VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
12919bd66eb0SPeter Brune       ierr = (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);CHKERRQ(ierr);
12929bd66eb0SPeter Brune     } else {
12936a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_func,   NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
12946a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_sol,    NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
12956a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
12966a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_func,     NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
12976a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_sol,      NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
12986a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_update,   NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
12996a388c36SPeter Brune     }
13009bd66eb0SPeter Brune   }
13016a388c36SPeter Brune   PetscFunctionReturn(0);
13026a388c36SPeter Brune }
13036a388c36SPeter Brune 
13046f263ca3SPeter Brune #undef __FUNCT__
13056f263ca3SPeter Brune #define __FUNCT__ "SNESLineSearchSetComputeNorms"
13066f263ca3SPeter Brune /*@
13076f263ca3SPeter Brune    SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search.
13086f263ca3SPeter Brune 
13096f263ca3SPeter Brune    Input Parameters:
131078bcb3b5SPeter Brune +  linesearch  - linesearch context
131178bcb3b5SPeter Brune -  flg  - indicates whether or not to compute norms
13126f263ca3SPeter Brune 
13136f263ca3SPeter Brune    Options Database Keys:
13148e557f58SPeter Brune .   -snes_linesearch_norms - turn norm computation on or off
13156f263ca3SPeter Brune 
13166f263ca3SPeter Brune    Notes:
13171a4f838cSPeter Brune    This is most relevant to the SNESLINESEARCHBASIC line search type.
13186f263ca3SPeter Brune 
13196f263ca3SPeter Brune    Level: intermediate
13206f263ca3SPeter Brune 
13211a4f838cSPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetNorms(), SNESLineSearchComputeNorms(), SNESLINESEARCHBASIC
13226f263ca3SPeter Brune @*/
13236f263ca3SPeter Brune PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg)
13246f263ca3SPeter Brune {
13256f263ca3SPeter Brune   PetscFunctionBegin;
13266f263ca3SPeter Brune   linesearch->norms = flg;
13276f263ca3SPeter Brune   PetscFunctionReturn(0);
13286f263ca3SPeter Brune }
13296f263ca3SPeter Brune 
13306a388c36SPeter Brune #undef __FUNCT__
1331f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetVecs"
1332f40b411bSPeter Brune /*@
1333f1c6b773SPeter Brune    SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context
1334f40b411bSPeter Brune 
1335f40b411bSPeter Brune    Input Parameters:
133678bcb3b5SPeter Brune .  linesearch - linesearch context
1337f40b411bSPeter Brune 
1338f40b411bSPeter Brune    Output Parameters:
13396232e825SPeter Brune +  X - Solution vector
13406232e825SPeter Brune .  F - Function vector
13416232e825SPeter Brune .  Y - Search direction vector
13426232e825SPeter Brune .  W - Solution work vector
13436232e825SPeter Brune -  G - Function work vector
13446232e825SPeter Brune 
13457bba9028SPeter Brune    Notes:
13467bba9028SPeter Brune    At the beginning of a line search application, X should contain a
13476232e825SPeter Brune    solution and the vector F the function computed at X.  At the end of the
13486232e825SPeter Brune    line search application, X should contain the new solution, and F the
13496232e825SPeter Brune    function evaluated at the new solution.
1350f40b411bSPeter Brune 
135178bcb3b5SPeter Brune    Level: advanced
1352f40b411bSPeter Brune 
1353f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1354f40b411bSPeter Brune @*/
1355bf388a1fSBarry Smith PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G)
1356bf388a1fSBarry Smith {
13576a388c36SPeter Brune   PetscFunctionBegin;
1358f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
13596a388c36SPeter Brune   if (X) {
13606a388c36SPeter Brune     PetscValidPointer(X, 2);
13616a388c36SPeter Brune     *X = linesearch->vec_sol;
13626a388c36SPeter Brune   }
13636a388c36SPeter Brune   if (F) {
13646a388c36SPeter Brune     PetscValidPointer(F, 3);
13656a388c36SPeter Brune     *F = linesearch->vec_func;
13666a388c36SPeter Brune   }
13676a388c36SPeter Brune   if (Y) {
13686a388c36SPeter Brune     PetscValidPointer(Y, 4);
13696a388c36SPeter Brune     *Y = linesearch->vec_update;
13706a388c36SPeter Brune   }
13716a388c36SPeter Brune   if (W) {
13726a388c36SPeter Brune     PetscValidPointer(W, 5);
13736a388c36SPeter Brune     *W = linesearch->vec_sol_new;
13746a388c36SPeter Brune   }
13756a388c36SPeter Brune   if (G) {
13766a388c36SPeter Brune     PetscValidPointer(G, 6);
13776a388c36SPeter Brune     *G = linesearch->vec_func_new;
13786a388c36SPeter Brune   }
13796a388c36SPeter Brune   PetscFunctionReturn(0);
13806a388c36SPeter Brune }
13816a388c36SPeter Brune 
13826a388c36SPeter Brune #undef __FUNCT__
1383f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetVecs"
1384f40b411bSPeter Brune /*@
1385f1c6b773SPeter Brune    SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context
1386f40b411bSPeter Brune 
1387f40b411bSPeter Brune    Input Parameters:
138878bcb3b5SPeter Brune +  linesearch - linesearch context
13896232e825SPeter Brune .  X - Solution vector
13906232e825SPeter Brune .  F - Function vector
13916232e825SPeter Brune .  Y - Search direction vector
13926232e825SPeter Brune .  W - Solution work vector
13936232e825SPeter Brune -  G - Function work vector
1394f40b411bSPeter Brune 
139578bcb3b5SPeter Brune    Level: advanced
1396f40b411bSPeter Brune 
1397f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms(), SNESLineSearchGetVecs()
1398f40b411bSPeter Brune @*/
1399bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G)
1400bf388a1fSBarry Smith {
14016a388c36SPeter Brune   PetscFunctionBegin;
1402f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
14036a388c36SPeter Brune   if (X) {
14046a388c36SPeter Brune     PetscValidHeaderSpecific(X,VEC_CLASSID,2);
14056a388c36SPeter Brune     linesearch->vec_sol = X;
14066a388c36SPeter Brune   }
14076a388c36SPeter Brune   if (F) {
14086a388c36SPeter Brune     PetscValidHeaderSpecific(F,VEC_CLASSID,3);
14096a388c36SPeter Brune     linesearch->vec_func = F;
14106a388c36SPeter Brune   }
14116a388c36SPeter Brune   if (Y) {
14126a388c36SPeter Brune     PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
14136a388c36SPeter Brune     linesearch->vec_update = Y;
14146a388c36SPeter Brune   }
14156a388c36SPeter Brune   if (W) {
14166a388c36SPeter Brune     PetscValidHeaderSpecific(W,VEC_CLASSID,5);
14176a388c36SPeter Brune     linesearch->vec_sol_new = W;
14186a388c36SPeter Brune   }
14196a388c36SPeter Brune   if (G) {
14206a388c36SPeter Brune     PetscValidHeaderSpecific(G,VEC_CLASSID,6);
14216a388c36SPeter Brune     linesearch->vec_func_new = G;
14226a388c36SPeter Brune   }
14236a388c36SPeter Brune   PetscFunctionReturn(0);
14246a388c36SPeter Brune }
14256a388c36SPeter Brune 
14266a388c36SPeter Brune #undef __FUNCT__
1427f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchAppendOptionsPrefix"
1428e7058c64SPeter Brune /*@C
1429f1c6b773SPeter Brune    SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1430e7058c64SPeter Brune    SNES options in the database.
1431e7058c64SPeter Brune 
1432cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
1433e7058c64SPeter Brune 
1434e7058c64SPeter Brune    Input Parameters:
1435e7058c64SPeter Brune +  snes - the SNES context
1436e7058c64SPeter Brune -  prefix - the prefix to prepend to all option names
1437e7058c64SPeter Brune 
1438e7058c64SPeter Brune    Notes:
1439e7058c64SPeter Brune    A hyphen (-) must NOT be given at the beginning of the prefix name.
1440e7058c64SPeter Brune    The first character of all runtime options is AUTOMATICALLY the hyphen.
1441e7058c64SPeter Brune 
1442e7058c64SPeter Brune    Level: advanced
1443e7058c64SPeter Brune 
1444f1c6b773SPeter Brune .keywords: SNESLineSearch, append, options, prefix, database
1445e7058c64SPeter Brune 
1446e7058c64SPeter Brune .seealso: SNESGetOptionsPrefix()
1447e7058c64SPeter Brune @*/
1448f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[])
1449e7058c64SPeter Brune {
1450e7058c64SPeter Brune   PetscErrorCode ierr;
1451e7058c64SPeter Brune 
1452e7058c64SPeter Brune   PetscFunctionBegin;
1453f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1454e7058c64SPeter Brune   ierr = PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1455e7058c64SPeter Brune   PetscFunctionReturn(0);
1456e7058c64SPeter Brune }
1457e7058c64SPeter Brune 
1458e7058c64SPeter Brune #undef __FUNCT__
1459f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetOptionsPrefix"
1460e7058c64SPeter Brune /*@C
1461f1c6b773SPeter Brune    SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all
1462f1c6b773SPeter Brune    SNESLineSearch options in the database.
1463e7058c64SPeter Brune 
1464e7058c64SPeter Brune    Not Collective
1465e7058c64SPeter Brune 
1466e7058c64SPeter Brune    Input Parameter:
1467cd7522eaSPeter Brune .  linesearch - the SNESLineSearch context
1468e7058c64SPeter Brune 
1469e7058c64SPeter Brune    Output Parameter:
1470e7058c64SPeter Brune .  prefix - pointer to the prefix string used
1471e7058c64SPeter Brune 
14728e557f58SPeter Brune    Notes:
14738e557f58SPeter Brune    On the fortran side, the user should pass in a string 'prefix' of
1474e7058c64SPeter Brune    sufficient length to hold the prefix.
1475e7058c64SPeter Brune 
1476e7058c64SPeter Brune    Level: advanced
1477e7058c64SPeter Brune 
1478f1c6b773SPeter Brune .keywords: SNESLineSearch, get, options, prefix, database
1479e7058c64SPeter Brune 
1480e7058c64SPeter Brune .seealso: SNESAppendOptionsPrefix()
1481e7058c64SPeter Brune @*/
1482f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[])
1483e7058c64SPeter Brune {
1484e7058c64SPeter Brune   PetscErrorCode ierr;
1485e7058c64SPeter Brune 
1486e7058c64SPeter Brune   PetscFunctionBegin;
1487f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1488e7058c64SPeter Brune   ierr = PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1489e7058c64SPeter Brune   PetscFunctionReturn(0);
1490e7058c64SPeter Brune }
1491bf7f4e0aSPeter Brune 
1492bf7f4e0aSPeter Brune #undef __FUNCT__
1493fa0ddf94SBarry Smith #define __FUNCT__ "SNESLineSearchSetWorkVecs"
14948d359177SBarry Smith /*@C
1495fa0ddf94SBarry Smith    SNESLineSearchSetWorkVecs - Gets work vectors for the line search.
1496f40b411bSPeter Brune 
1497f40b411bSPeter Brune    Input Parameter:
1498f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
1499f40b411bSPeter Brune -  nwork - the number of work vectors
1500f40b411bSPeter Brune 
1501f40b411bSPeter Brune    Level: developer
1502f40b411bSPeter Brune 
1503fa0ddf94SBarry Smith    Developers Note: This is PETSC_EXTERN because it may be used by user written plugin SNES implementations
1504cd7522eaSPeter Brune 
1505f1c6b773SPeter Brune .keywords: SNESLineSearch, work, vector
1506f40b411bSPeter Brune 
1507fa0ddf94SBarry Smith .seealso: SNESSetWorkVecs()
1508f40b411bSPeter Brune @*/
1509fa0ddf94SBarry Smith PetscErrorCode  SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork)
1510bf7f4e0aSPeter Brune {
1511bf7f4e0aSPeter Brune   PetscErrorCode ierr;
1512bf388a1fSBarry Smith 
1513bf7f4e0aSPeter Brune   PetscFunctionBegin;
1514bf7f4e0aSPeter Brune   if (linesearch->vec_sol) {
1515bf7f4e0aSPeter Brune     ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr);
15168d359177SBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
1517bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1518bf7f4e0aSPeter Brune }
1519bf7f4e0aSPeter Brune 
1520bf7f4e0aSPeter Brune #undef __FUNCT__
1521f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetSuccess"
1522f40b411bSPeter Brune /*@
1523f1c6b773SPeter Brune    SNESLineSearchGetSuccess - Gets the success/failure status of the last line search application
1524f40b411bSPeter Brune 
1525f40b411bSPeter Brune    Input Parameters:
152678bcb3b5SPeter Brune .  linesearch - linesearch context
1527f40b411bSPeter Brune 
1528f40b411bSPeter Brune    Output Parameters:
15298e557f58SPeter Brune .  success - The success or failure status
1530f40b411bSPeter Brune 
1531cd7522eaSPeter Brune    Notes:
1532c98378a5SBarry Smith    This is typically called after SNESLineSearchApply() in order to determine if the line-search failed
1533cd7522eaSPeter Brune    (and set the SNES convergence accordingly).
1534cd7522eaSPeter Brune 
1535f40b411bSPeter Brune    Level: intermediate
1536f40b411bSPeter Brune 
1537f1c6b773SPeter Brune .seealso: SNESLineSearchSetSuccess()
1538f40b411bSPeter Brune @*/
1539f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetSuccess(SNESLineSearch linesearch, PetscBool *success)
1540bf7f4e0aSPeter Brune {
1541bf7f4e0aSPeter Brune   PetscFunctionBegin;
1542f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
15436a388c36SPeter Brune   PetscValidPointer(success, 2);
1544f5af7f23SKarl Rupp   if (success) *success = linesearch->success;
1545bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1546bf7f4e0aSPeter Brune }
1547bf7f4e0aSPeter Brune 
1548bf7f4e0aSPeter Brune #undef __FUNCT__
1549f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetSuccess"
1550f40b411bSPeter Brune /*@
1551f1c6b773SPeter Brune    SNESLineSearchSetSuccess - Sets the success/failure status of the last line search application
1552f40b411bSPeter Brune 
1553f40b411bSPeter Brune    Input Parameters:
155478bcb3b5SPeter Brune +  linesearch - linesearch context
15558e557f58SPeter Brune -  success - The success or failure status
1556f40b411bSPeter Brune 
1557cd7522eaSPeter Brune    Notes:
1558cd7522eaSPeter Brune    This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set
1559cd7522eaSPeter Brune    the success or failure of the line search method.
1560cd7522eaSPeter Brune 
156178bcb3b5SPeter Brune    Level: developer
1562f40b411bSPeter Brune 
1563f1c6b773SPeter Brune .seealso: SNESLineSearchGetSuccess()
1564f40b411bSPeter Brune @*/
1565f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetSuccess(SNESLineSearch linesearch, PetscBool success)
15666a388c36SPeter Brune {
15676a388c36SPeter Brune   PetscFunctionBegin;
15685fd66863SKarl Rupp   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
15696a388c36SPeter Brune   linesearch->success = success;
15706a388c36SPeter Brune   PetscFunctionReturn(0);
15716a388c36SPeter Brune }
15726a388c36SPeter Brune 
15736a388c36SPeter Brune #undef __FUNCT__
1574f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetVIFunctions"
15759bd66eb0SPeter Brune /*@C
1576f1c6b773SPeter Brune    SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.
15779bd66eb0SPeter Brune 
15789bd66eb0SPeter Brune    Input Parameters:
15799bd66eb0SPeter Brune +  snes - nonlinear context obtained from SNESCreate()
15809bd66eb0SPeter Brune .  projectfunc - function for projecting the function to the bounds
15819bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
15829bd66eb0SPeter Brune 
15839bd66eb0SPeter Brune    Logically Collective on SNES
15849bd66eb0SPeter Brune 
15859bd66eb0SPeter Brune    Calling sequence of projectfunc:
15869bd66eb0SPeter Brune .vb
15879bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X)
15889bd66eb0SPeter Brune .ve
15899bd66eb0SPeter Brune 
15909bd66eb0SPeter Brune     Input parameters for projectfunc:
15919bd66eb0SPeter Brune +   snes - nonlinear context
15929bd66eb0SPeter Brune -   X - current solution
15939bd66eb0SPeter Brune 
1594cd7522eaSPeter Brune     Output parameters for projectfunc:
15959bd66eb0SPeter Brune .   X - Projected solution
15969bd66eb0SPeter Brune 
15979bd66eb0SPeter Brune    Calling sequence of normfunc:
15989bd66eb0SPeter Brune .vb
15999bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm)
16009bd66eb0SPeter Brune .ve
16019bd66eb0SPeter Brune 
1602cd7522eaSPeter Brune     Input parameters for normfunc:
16039bd66eb0SPeter Brune +   snes - nonlinear context
16049bd66eb0SPeter Brune .   X - current solution
16059bd66eb0SPeter Brune -   F - current residual
16069bd66eb0SPeter Brune 
1607cd7522eaSPeter Brune     Output parameters for normfunc:
16089bd66eb0SPeter Brune .   fnorm - VI-specific norm of the function
16099bd66eb0SPeter Brune 
1610cd7522eaSPeter Brune     Notes:
1611cd7522eaSPeter Brune     The VI solvers require projection of the solution to the feasible set.  projectfunc should implement this.
1612cd7522eaSPeter Brune 
1613cd7522eaSPeter Brune     The VI solvers require special evaluation of the function norm such that the norm is only calculated
1614cd7522eaSPeter Brune     on the inactive set.  This should be implemented by normfunc.
16159bd66eb0SPeter Brune 
16169bd66eb0SPeter Brune     Level: developer
16179bd66eb0SPeter Brune 
16189bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, set, line search
16199bd66eb0SPeter Brune 
1620f1c6b773SPeter Brune .seealso: SNESLineSearchGetVIFunctions(), SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck()
16219bd66eb0SPeter Brune @*/
1622f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc)
16239bd66eb0SPeter Brune {
16249bd66eb0SPeter Brune   PetscFunctionBegin;
1625f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
16269bd66eb0SPeter Brune   if (projectfunc) linesearch->ops->viproject = projectfunc;
16279bd66eb0SPeter Brune   if (normfunc) linesearch->ops->vinorm = normfunc;
16289bd66eb0SPeter Brune   PetscFunctionReturn(0);
16299bd66eb0SPeter Brune }
16309bd66eb0SPeter Brune 
16319bd66eb0SPeter Brune #undef __FUNCT__
1632f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetVIFunctions"
16339bd66eb0SPeter Brune /*@C
1634f1c6b773SPeter Brune    SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.
16359bd66eb0SPeter Brune 
16369bd66eb0SPeter Brune    Input Parameters:
16379bd66eb0SPeter Brune .  snes - nonlinear context obtained from SNESCreate()
16389bd66eb0SPeter Brune 
16399bd66eb0SPeter Brune    Output Parameters:
16409bd66eb0SPeter Brune +  projectfunc - function for projecting the function to the bounds
16419bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
16429bd66eb0SPeter Brune 
16439bd66eb0SPeter Brune    Logically Collective on SNES
16449bd66eb0SPeter Brune 
16459bd66eb0SPeter Brune     Level: developer
16469bd66eb0SPeter Brune 
16479bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, get, line search
16489bd66eb0SPeter Brune 
1649f1c6b773SPeter Brune .seealso: SNESLineSearchSetVIFunctions(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
16509bd66eb0SPeter Brune @*/
1651f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc)
16529bd66eb0SPeter Brune {
16539bd66eb0SPeter Brune   PetscFunctionBegin;
16549bd66eb0SPeter Brune   if (projectfunc) *projectfunc = linesearch->ops->viproject;
16559bd66eb0SPeter Brune   if (normfunc) *normfunc = linesearch->ops->vinorm;
16569bd66eb0SPeter Brune   PetscFunctionReturn(0);
16579bd66eb0SPeter Brune }
16589bd66eb0SPeter Brune 
16599bd66eb0SPeter Brune #undef __FUNCT__
1660f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchRegister"
1661bf7f4e0aSPeter Brune /*@C
16621c84c290SBarry Smith   SNESLineSearchRegister - See SNESLineSearchRegister()
1663bf7f4e0aSPeter Brune 
1664bf7f4e0aSPeter Brune   Level: advanced
1665bf7f4e0aSPeter Brune @*/
1666bdf89e91SBarry Smith PetscErrorCode  SNESLineSearchRegister(const char sname[],PetscErrorCode (*function)(SNESLineSearch))
1667bf7f4e0aSPeter Brune {
1668bf7f4e0aSPeter Brune   PetscErrorCode ierr;
1669bf7f4e0aSPeter Brune 
1670bf7f4e0aSPeter Brune   PetscFunctionBegin;
1671a240a19fSJed Brown   ierr = PetscFunctionListAdd(&SNESLineSearchList,sname,function);CHKERRQ(ierr);
1672bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1673bf7f4e0aSPeter Brune }
1674