xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision 7601faf09de2bc280184d9c73acff348d5eb2a25)
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);
360298fd71SBarry Smith   *outlinesearch = NULL;
37f5af7f23SKarl Rupp 
3867c2884eSBarry Smith   ierr = PetscHeaderCreate(linesearch,_p_LineSearch,struct _LineSearchOps,SNESLINESEARCH_CLASSID, "SNESLineSearch","Linesearch","SNESLineSearch",comm,SNESLineSearchDestroy,SNESLineSearchView);CHKERRQ(ierr);
39bf7f4e0aSPeter Brune 
400298fd71SBarry Smith   linesearch->ops->precheck  = NULL;
410298fd71SBarry Smith   linesearch->ops->postcheck = NULL;
42bf7f4e0aSPeter Brune 
430298fd71SBarry Smith   linesearch->vec_sol_new  = NULL;
440298fd71SBarry Smith   linesearch->vec_func_new = NULL;
450298fd71SBarry Smith   linesearch->vec_sol      = NULL;
460298fd71SBarry Smith   linesearch->vec_func     = NULL;
470298fd71SBarry Smith   linesearch->vec_update   = NULL;
489bd66eb0SPeter Brune 
49bf7f4e0aSPeter Brune   linesearch->lambda       = 1.0;
50bf7f4e0aSPeter Brune   linesearch->fnorm        = 1.0;
51bf7f4e0aSPeter Brune   linesearch->ynorm        = 1.0;
52bf7f4e0aSPeter Brune   linesearch->xnorm        = 1.0;
53bf7f4e0aSPeter Brune   linesearch->success      = PETSC_TRUE;
54bf7f4e0aSPeter Brune   linesearch->norms        = PETSC_TRUE;
55bf7f4e0aSPeter Brune   linesearch->keeplambda   = PETSC_FALSE;
56bf7f4e0aSPeter Brune   linesearch->damping      = 1.0;
57bf7f4e0aSPeter Brune   linesearch->maxstep      = 1e8;
58bf7f4e0aSPeter Brune   linesearch->steptol      = 1e-12;
59516fe3c3SPeter Brune   linesearch->rtol         = 1e-8;
60516fe3c3SPeter Brune   linesearch->atol         = 1e-15;
61516fe3c3SPeter Brune   linesearch->ltol         = 1e-8;
620298fd71SBarry Smith   linesearch->precheckctx  = NULL;
630298fd71SBarry Smith   linesearch->postcheckctx = NULL;
6459405d5eSPeter Brune   linesearch->max_its      = 1;
65bf7f4e0aSPeter Brune   linesearch->setupcalled  = PETSC_FALSE;
66bf7f4e0aSPeter Brune   *outlinesearch           = linesearch;
67bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
68bf7f4e0aSPeter Brune }
69bf7f4e0aSPeter Brune 
70bf7f4e0aSPeter Brune #undef __FUNCT__
71f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetUp"
72f40b411bSPeter Brune /*@
7378bcb3b5SPeter Brune    SNESLineSearchSetUp - Prepares the line search for being applied by allocating
7478bcb3b5SPeter Brune    any required vectors.
75f40b411bSPeter Brune 
76cd7522eaSPeter Brune    Collective on SNESLineSearch
77f40b411bSPeter Brune 
78f40b411bSPeter Brune    Input Parameters:
79f40b411bSPeter Brune .  linesearch - The LineSearch instance.
80f40b411bSPeter Brune 
81cd7522eaSPeter Brune    Notes:
82cd7522eaSPeter Brune    For most cases, this needn't be called outside of SNESLineSearchApply().
83cd7522eaSPeter Brune    The only current case where this is called outside of this is for the VI
8478bcb3b5SPeter Brune    solvers, which modify the solution and work vectors before the first call
85cd7522eaSPeter Brune    of SNESLineSearchApply, requiring the SNESLineSearch work vectors to be
86cd7522eaSPeter Brune    allocated upfront.
87cd7522eaSPeter Brune 
8878bcb3b5SPeter Brune    Level: advanced
89f40b411bSPeter Brune 
90f1c6b773SPeter Brune .keywords: SNESLineSearch, SetUp
91f40b411bSPeter Brune 
92f1c6b773SPeter Brune .seealso: SNESLineSearchReset()
93f40b411bSPeter Brune @*/
94f40b411bSPeter Brune 
95bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch)
96bf388a1fSBarry Smith {
97bf7f4e0aSPeter Brune   PetscErrorCode ierr;
98bf388a1fSBarry Smith 
99bf7f4e0aSPeter Brune   PetscFunctionBegin;
100bf7f4e0aSPeter Brune   if (!((PetscObject)linesearch)->type_name) {
1011a4f838cSPeter Brune     ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr);
102bf7f4e0aSPeter Brune   }
103bf7f4e0aSPeter Brune   if (!linesearch->setupcalled) {
1049bd66eb0SPeter Brune     if (!linesearch->vec_sol_new) {
105bf7f4e0aSPeter Brune       ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new);CHKERRQ(ierr);
1069bd66eb0SPeter Brune     }
1079bd66eb0SPeter Brune     if (!linesearch->vec_func_new) {
1089bd66eb0SPeter Brune       ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new);CHKERRQ(ierr);
1099bd66eb0SPeter Brune     }
110bf7f4e0aSPeter Brune     if (linesearch->ops->setup) {
111bf7f4e0aSPeter Brune       ierr = (*linesearch->ops->setup)(linesearch);CHKERRQ(ierr);
112bf7f4e0aSPeter Brune     }
113bf7f4e0aSPeter Brune     linesearch->lambda      = linesearch->damping;
114bf7f4e0aSPeter Brune     linesearch->setupcalled = PETSC_TRUE;
115bf7f4e0aSPeter Brune   }
116bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
117bf7f4e0aSPeter Brune }
118bf7f4e0aSPeter Brune 
119bf7f4e0aSPeter Brune #undef __FUNCT__
120f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchReset"
121f40b411bSPeter Brune 
122f40b411bSPeter Brune /*@
12378bcb3b5SPeter Brune    SNESLineSearchReset - Undoes the SetUp and deletes any Vecs or Mats allocated by the line search.
124f40b411bSPeter Brune 
125f1c6b773SPeter Brune    Collective on SNESLineSearch
126f40b411bSPeter Brune 
127f40b411bSPeter Brune    Input Parameters:
128f40b411bSPeter Brune .  linesearch - The LineSearch instance.
129f40b411bSPeter Brune 
13078bcb3b5SPeter Brune    Level: intermediate
131f40b411bSPeter Brune 
132cd7522eaSPeter Brune .keywords: SNESLineSearch, Reset
133f40b411bSPeter Brune 
134f1c6b773SPeter Brune .seealso: SNESLineSearchSetUp()
135f40b411bSPeter Brune @*/
136f40b411bSPeter Brune 
137bf388a1fSBarry Smith PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch)
138bf388a1fSBarry Smith {
139bf7f4e0aSPeter Brune   PetscErrorCode ierr;
140bf388a1fSBarry Smith 
141bf7f4e0aSPeter Brune   PetscFunctionBegin;
142f5af7f23SKarl Rupp   if (linesearch->ops->reset) (*linesearch->ops->reset)(linesearch);
143f5af7f23SKarl Rupp 
144bf7f4e0aSPeter Brune   ierr = VecDestroy(&linesearch->vec_sol_new);CHKERRQ(ierr);
145bf7f4e0aSPeter Brune   ierr = VecDestroy(&linesearch->vec_func_new);CHKERRQ(ierr);
146bf7f4e0aSPeter Brune 
147bf7f4e0aSPeter Brune   ierr = VecDestroyVecs(linesearch->nwork, &linesearch->work);CHKERRQ(ierr);
148f5af7f23SKarl Rupp 
149bf7f4e0aSPeter Brune   linesearch->nwork       = 0;
150bf7f4e0aSPeter Brune   linesearch->setupcalled = PETSC_FALSE;
151bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
152bf7f4e0aSPeter Brune }
153bf7f4e0aSPeter Brune 
1546b2b7091SBarry Smith /*MC
1556b2b7091SBarry Smith     SNESLineSearchPreCheckFunction - functional form passed to check before line search is called
1566b2b7091SBarry Smith 
1576b2b7091SBarry Smith      Synopsis:
1586b2b7091SBarry Smith      #include "petscsnes.h"
1596b2b7091SBarry Smith      SNESLineSearchPreCheckFunction(SNESLineSearch snes,Vec x,Vec y, PetscBool *changed);
1606b2b7091SBarry Smith 
1616b2b7091SBarry Smith        Input Parameters:
1626b2b7091SBarry Smith +      x - solution vector
1636b2b7091SBarry Smith .      y - search direction vector
1646b2b7091SBarry Smith -      changed - flag to indicate the precheck changed x or y.
1656b2b7091SBarry Smith 
166878cb397SSatish Balay    Level: advanced
167878cb397SSatish Balay 
1686b2b7091SBarry Smith .seealso:   SNESLineSearchSetPreCheck(), SNESLineSearchSetPostCheck()
1696b2b7091SBarry Smith M*/
1706b2b7091SBarry Smith 
17186d74e61SPeter Brune #undef __FUNCT__
172f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetPreCheck"
17386d74e61SPeter Brune /*@C
174f1c6b773SPeter Brune    SNESLineSearchSetPreCheck - Sets a pre-check function for the line search routine.
17586d74e61SPeter Brune 
176f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
17786d74e61SPeter Brune 
17886d74e61SPeter Brune    Input Parameters:
179f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
1806b2b7091SBarry Smith .  SNESLineSearchPreCheckFunction - [optional] function evaluation routine
18186d74e61SPeter Brune -  ctx        - [optional] user-defined context for private data for the
1820298fd71SBarry Smith                 function evaluation routine (may be NULL)
18386d74e61SPeter Brune 
18486d74e61SPeter Brune 
18586d74e61SPeter Brune    Level: intermediate
18686d74e61SPeter Brune 
18786d74e61SPeter Brune .keywords: set, linesearch, pre-check
18886d74e61SPeter Brune 
189f1c6b773SPeter Brune .seealso: SNESLineSearchSetPostCheck()
19086d74e61SPeter Brune @*/
1916b2b7091SBarry Smith PetscErrorCode  SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*SNESLineSearchPreCheckFunction)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void *ctx)
19286d74e61SPeter Brune {
1939bd66eb0SPeter Brune   PetscFunctionBegin;
194f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1956b2b7091SBarry Smith   if (SNESLineSearchPreCheckFunction) linesearch->ops->precheck = SNESLineSearchPreCheckFunction;
19686d74e61SPeter Brune   if (ctx) linesearch->precheckctx = ctx;
19786d74e61SPeter Brune   PetscFunctionReturn(0);
19886d74e61SPeter Brune }
19986d74e61SPeter Brune 
20086d74e61SPeter Brune #undef __FUNCT__
201f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetPreCheck"
20286d74e61SPeter Brune /*@C
203cd7522eaSPeter Brune    SNESLineSearchSetPreCheck - Gets the pre-check function for the line search routine.
20486d74e61SPeter Brune 
20586d74e61SPeter Brune    Input Parameters:
206f1c6b773SPeter Brune .  linesearch - the SNESLineSearch context
20786d74e61SPeter Brune 
20886d74e61SPeter Brune    Output Parameters:
20986d74e61SPeter Brune +  func       - [optional] function evaluation routine
21086d74e61SPeter Brune -  ctx        - [optional] user-defined context for private data for the
2110298fd71SBarry Smith                 function evaluation routine (may be NULL)
21286d74e61SPeter Brune 
21386d74e61SPeter Brune    Level: intermediate
21486d74e61SPeter Brune 
21586d74e61SPeter Brune .keywords: get, linesearch, pre-check
21686d74e61SPeter Brune 
217f1c6b773SPeter Brune .seealso: SNESLineSearchGetPostCheck(), SNESLineSearchSetPreCheck()
21886d74e61SPeter Brune @*/
2196b2b7091SBarry Smith PetscErrorCode  SNESLineSearchGetPreCheck(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) *SNESLineSearchPreCheckFunction = linesearch->ops->precheck;
22486d74e61SPeter Brune   if (ctx) *ctx = linesearch->precheckctx;
22586d74e61SPeter Brune   PetscFunctionReturn(0);
22686d74e61SPeter Brune }
22786d74e61SPeter Brune 
2286b2b7091SBarry Smith /*MC
2296b2b7091SBarry Smith     SNESLineSearchPostheckFunction - functional form that is called after line search is complete
2306b2b7091SBarry Smith 
2316b2b7091SBarry Smith      Synopsis:
2326b2b7091SBarry Smith      #include "petscsnes.h"
2336b2b7091SBarry Smith      SNESLineSearchPostheckFunction(SNESLineSearch linesearch,Vec x,Vec y,  Vec w, *changed_y, PetscBool *changed_w);
2346b2b7091SBarry Smith 
2356b2b7091SBarry Smith      Input Parameters:
2366b2b7091SBarry Smith +      x - old solution vector
2376b2b7091SBarry Smith .      y - search direction vector
2386b2b7091SBarry Smith .      w - new solution vector
2396b2b7091SBarry Smith .      changed_y - indicates that the line search changed y
2406b2b7091SBarry Smith -      changed_w - indicates that the line search changed w
2416b2b7091SBarry Smith 
242878cb397SSatish Balay    Level: advanced
2436b2b7091SBarry Smith 
2446b2b7091SBarry Smith .seealso:   SNESLineSearchSetPreCheck(), SNESLineSearchSetPostCheck()
2456b2b7091SBarry Smith M*/
2466b2b7091SBarry Smith 
24786d74e61SPeter Brune #undef __FUNCT__
248f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetPostCheck"
24986d74e61SPeter Brune /*@C
250f1c6b773SPeter Brune    SNESLineSearchSetPostCheck - Sets a post-check function for the line search routine.
25186d74e61SPeter Brune 
252f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
25386d74e61SPeter Brune 
25486d74e61SPeter Brune    Input Parameters:
255f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
2566b2b7091SBarry Smith .  SNESLineSearchPostCheckFunction - [optional] function evaluation routine
25786d74e61SPeter Brune -  ctx        - [optional] user-defined context for private data for the
2580298fd71SBarry Smith                 function evaluation routine (may be NULL)
25986d74e61SPeter Brune 
26086d74e61SPeter Brune    Level: intermediate
26186d74e61SPeter Brune 
26286d74e61SPeter Brune .keywords: set, linesearch, post-check
26386d74e61SPeter Brune 
264f1c6b773SPeter Brune .seealso: SNESLineSearchSetPreCheck()
26586d74e61SPeter Brune @*/
2666b2b7091SBarry Smith PetscErrorCode  SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*SNESLineSearchPostCheckFunction)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void *ctx)
26786d74e61SPeter Brune {
26886d74e61SPeter Brune   PetscFunctionBegin;
269f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
2706b2b7091SBarry Smith   if (SNESLineSearchPostCheckFunction) linesearch->ops->postcheck = SNESLineSearchPostCheckFunction;
27186d74e61SPeter Brune   if (ctx) linesearch->postcheckctx = ctx;
27286d74e61SPeter Brune   PetscFunctionReturn(0);
27386d74e61SPeter Brune }
27486d74e61SPeter Brune 
27586d74e61SPeter Brune #undef __FUNCT__
276f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetPostCheck"
27786d74e61SPeter Brune /*@C
278f1c6b773SPeter Brune    SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine.
27986d74e61SPeter Brune 
28086d74e61SPeter Brune    Input Parameters:
281f1c6b773SPeter Brune .  linesearch - the SNESLineSearch context
28286d74e61SPeter Brune 
28386d74e61SPeter Brune    Output Parameters:
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: get, linesearch, post-check
29186d74e61SPeter Brune 
292f1c6b773SPeter Brune .seealso: SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck()
29386d74e61SPeter Brune @*/
2946b2b7091SBarry Smith PetscErrorCode  SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**SNESLineSearchPostCheckFunction)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void **ctx)
29586d74e61SPeter Brune {
2969bd66eb0SPeter Brune   PetscFunctionBegin;
297f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
2986b2b7091SBarry Smith   if (SNESLineSearchPostCheckFunction) *SNESLineSearchPostCheckFunction = linesearch->ops->postcheck;
29986d74e61SPeter Brune   if (ctx) *ctx = linesearch->postcheckctx;
30086d74e61SPeter Brune   PetscFunctionReturn(0);
30186d74e61SPeter Brune }
30286d74e61SPeter Brune 
303bf7f4e0aSPeter Brune #undef __FUNCT__
304f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPreCheck"
305f40b411bSPeter Brune /*@
306f1c6b773SPeter Brune    SNESLineSearchPreCheck - Prepares the line search for being applied.
307f40b411bSPeter Brune 
308cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
309f40b411bSPeter Brune 
310f40b411bSPeter Brune    Input Parameters:
3117b1df9c1SPeter Brune +  linesearch - The linesearch instance.
3127b1df9c1SPeter Brune .  X - The current solution
3137b1df9c1SPeter Brune -  Y - The step direction
314f40b411bSPeter Brune 
315f40b411bSPeter Brune    Output Parameters:
3168e557f58SPeter Brune .  changed - Indicator that the precheck routine has changed anything
317f40b411bSPeter Brune 
318f40b411bSPeter Brune    Level: Beginner
319f40b411bSPeter Brune 
320f1c6b773SPeter Brune .keywords: SNESLineSearch, Create
321f40b411bSPeter Brune 
322f1c6b773SPeter Brune .seealso: SNESLineSearchPostCheck()
323f40b411bSPeter Brune @*/
3247b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed)
325bf7f4e0aSPeter Brune {
326bf7f4e0aSPeter Brune   PetscErrorCode ierr;
3275fd66863SKarl Rupp 
328bf7f4e0aSPeter Brune   PetscFunctionBegin;
329bf7f4e0aSPeter Brune   *changed = PETSC_FALSE;
3306b2b7091SBarry Smith   if (linesearch->ops->precheck) {
3316b2b7091SBarry Smith     ierr = (*linesearch->ops->precheck)(linesearch, X, Y, changed, linesearch->precheckctx);CHKERRQ(ierr);
33238bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed,4);
333bf7f4e0aSPeter Brune   }
334bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
335bf7f4e0aSPeter Brune }
336bf7f4e0aSPeter Brune 
337bf7f4e0aSPeter Brune #undef __FUNCT__
338f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPostCheck"
339f40b411bSPeter Brune /*@
340f1c6b773SPeter Brune    SNESLineSearchPostCheck - Prepares the line search for being applied.
341f40b411bSPeter Brune 
342cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
343f40b411bSPeter Brune 
344f40b411bSPeter Brune    Input Parameters:
3457b1df9c1SPeter Brune +  linesearch - The linesearch context
3467b1df9c1SPeter Brune .  X - The last solution
3477b1df9c1SPeter Brune .  Y - The step direction
3487b1df9c1SPeter Brune -  W - The updated solution, W = X + lambda*Y for some lambda
349f40b411bSPeter Brune 
350f40b411bSPeter Brune    Output Parameters:
35178bcb3b5SPeter Brune +  changed_Y - Indicator if the direction Y has been changed.
35278bcb3b5SPeter Brune -  changed_W - Indicator if the new candidate solution W has been changed.
353f40b411bSPeter Brune 
354f40b411bSPeter Brune    Level: Intermediate
355f40b411bSPeter Brune 
356f1c6b773SPeter Brune .keywords: SNESLineSearch, Create
357f40b411bSPeter Brune 
358f1c6b773SPeter Brune .seealso: SNESLineSearchPreCheck()
359f40b411bSPeter Brune @*/
3607b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch,Vec X,Vec Y,Vec W,PetscBool *changed_Y,PetscBool *changed_W)
361bf7f4e0aSPeter Brune {
362bf7f4e0aSPeter Brune   PetscErrorCode ierr;
363bf388a1fSBarry Smith 
364bf7f4e0aSPeter Brune   PetscFunctionBegin;
365bf7f4e0aSPeter Brune   *changed_Y = PETSC_FALSE;
366bf7f4e0aSPeter Brune   *changed_W = PETSC_FALSE;
3676b2b7091SBarry Smith   if (linesearch->ops->postcheck) {
3686b2b7091SBarry Smith     ierr = (*linesearch->ops->postcheck)(linesearch,X,Y,W,changed_Y,changed_W,linesearch->postcheckctx);CHKERRQ(ierr);
36938bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed_Y,5);
37038bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed_W,6);
37186d74e61SPeter Brune   }
37286d74e61SPeter Brune   PetscFunctionReturn(0);
37386d74e61SPeter Brune }
37486d74e61SPeter Brune 
37586d74e61SPeter Brune #undef __FUNCT__
376f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPreCheckPicard"
37786d74e61SPeter Brune /*@C
37886d74e61SPeter Brune    SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration
37986d74e61SPeter Brune 
380cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
38186d74e61SPeter Brune 
38286d74e61SPeter Brune    Input Arguments:
38386d74e61SPeter Brune +  linesearch - linesearch context
38486d74e61SPeter Brune .  X - base state for this step
38586d74e61SPeter Brune .  Y - initial correction
38686d74e61SPeter Brune 
38786d74e61SPeter Brune    Output Arguments:
38886d74e61SPeter Brune +  Y - correction, possibly modified
38986d74e61SPeter Brune -  changed - flag indicating that Y was modified
39086d74e61SPeter Brune 
39186d74e61SPeter Brune    Options Database Key:
392cd7522eaSPeter Brune +  -snes_linesearch_precheck_picard - activate this routine
393cd7522eaSPeter Brune -  -snes_linesearch_precheck_picard_angle - angle
39486d74e61SPeter Brune 
39586d74e61SPeter Brune    Level: advanced
39686d74e61SPeter Brune 
39786d74e61SPeter Brune    Notes:
39886d74e61SPeter Brune    This function should be passed to SNESLineSearchSetPreCheck()
39986d74e61SPeter Brune 
40086d74e61SPeter Brune    The justification for this method involves the linear convergence of a Picard iteration
40186d74e61SPeter Brune    so the Picard linearization should be provided in place of the "Jacobian". This correction
40286d74e61SPeter Brune    is generally not useful when using a Newton linearization.
40386d74e61SPeter Brune 
40486d74e61SPeter Brune    Reference:
40586d74e61SPeter Brune    Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology.
40686d74e61SPeter Brune 
40786d74e61SPeter Brune .seealso: SNESLineSearchSetPreCheck()
40886d74e61SPeter Brune @*/
409f1c6b773SPeter Brune PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed,void *ctx)
41086d74e61SPeter Brune {
41186d74e61SPeter Brune   PetscErrorCode ierr;
41286d74e61SPeter Brune   PetscReal      angle = *(PetscReal*)linesearch->precheckctx;
41386d74e61SPeter Brune   Vec            Ylast;
41486d74e61SPeter Brune   PetscScalar    dot;
41586d74e61SPeter Brune   PetscInt       iter;
41686d74e61SPeter Brune   PetscReal      ynorm,ylastnorm,theta,angle_radians;
41786d74e61SPeter Brune   SNES           snes;
41886d74e61SPeter Brune 
41986d74e61SPeter Brune   PetscFunctionBegin;
420f1c6b773SPeter Brune   ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
42186d74e61SPeter Brune   ierr = PetscObjectQuery((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject*)&Ylast);CHKERRQ(ierr);
42286d74e61SPeter Brune   if (!Ylast) {
42386d74e61SPeter Brune     ierr = VecDuplicate(Y,&Ylast);CHKERRQ(ierr);
42486d74e61SPeter Brune     ierr = PetscObjectCompose((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject)Ylast);CHKERRQ(ierr);
42586d74e61SPeter Brune     ierr = PetscObjectDereference((PetscObject)Ylast);CHKERRQ(ierr);
42686d74e61SPeter Brune   }
42786d74e61SPeter Brune   ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr);
42886d74e61SPeter Brune   if (iter < 2) {
42986d74e61SPeter Brune     ierr     = VecCopy(Y,Ylast);CHKERRQ(ierr);
43086d74e61SPeter Brune     *changed = PETSC_FALSE;
43186d74e61SPeter Brune     PetscFunctionReturn(0);
43286d74e61SPeter Brune   }
43386d74e61SPeter Brune 
43486d74e61SPeter Brune   ierr = VecDot(Y,Ylast,&dot);CHKERRQ(ierr);
43586d74e61SPeter Brune   ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr);
43686d74e61SPeter Brune   ierr = VecNorm(Ylast,NORM_2,&ylastnorm);CHKERRQ(ierr);
43786d74e61SPeter Brune   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
43886d74e61SPeter Brune   theta         = acos((double)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0));
43986d74e61SPeter Brune   angle_radians = angle * PETSC_PI / 180.;
44086d74e61SPeter Brune   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
44186d74e61SPeter Brune     /* Modify the step Y */
44286d74e61SPeter Brune     PetscReal alpha,ydiffnorm;
44386d74e61SPeter Brune     ierr  = VecAXPY(Ylast,-1.0,Y);CHKERRQ(ierr);
44486d74e61SPeter Brune     ierr  = VecNorm(Ylast,NORM_2,&ydiffnorm);CHKERRQ(ierr);
44586d74e61SPeter Brune     alpha = ylastnorm / ydiffnorm;
44686d74e61SPeter Brune     ierr  = VecCopy(Y,Ylast);CHKERRQ(ierr);
44786d74e61SPeter Brune     ierr  = VecScale(Y,alpha);CHKERRQ(ierr);
448c69d1a72SBarry 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);
44986d74e61SPeter Brune   } else {
450c69d1a72SBarry 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);
45186d74e61SPeter Brune     ierr     = VecCopy(Y,Ylast);CHKERRQ(ierr);
45286d74e61SPeter Brune     *changed = PETSC_FALSE;
453bf7f4e0aSPeter Brune   }
454bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
455bf7f4e0aSPeter Brune }
456bf7f4e0aSPeter Brune 
457bf7f4e0aSPeter Brune #undef __FUNCT__
458f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchApply"
459f40b411bSPeter Brune /*@
460cd7522eaSPeter Brune    SNESLineSearchApply - Computes the line-search update.
461f40b411bSPeter Brune 
462f1c6b773SPeter Brune    Collective on SNESLineSearch
463f40b411bSPeter Brune 
464f40b411bSPeter Brune    Input Parameters:
4658e557f58SPeter Brune +  linesearch - The linesearch context
4668e557f58SPeter Brune .  X - The current solution
4678e557f58SPeter Brune .  F - The current function
4688e557f58SPeter Brune .  fnorm - The current norm
4698e557f58SPeter Brune -  Y - The search direction
470f40b411bSPeter Brune 
471f40b411bSPeter Brune    Output Parameters:
4728e557f58SPeter Brune +  X - The new solution
4738e557f58SPeter Brune .  F - The new function
4748e557f58SPeter Brune -  fnorm - The new function norm
475f40b411bSPeter Brune 
476cd7522eaSPeter Brune    Options Database Keys:
4773c7d6663SPeter Brune + -snes_linesearch_type - basic, bt, l2, cp, shell
478cd7522eaSPeter Brune . -snes_linesearch_monitor - Print progress of line searches
4793c7d6663SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter
480cd7522eaSPeter Brune . -snes_linesearch_norms   - Turn on/off the linesearch norms
4813c7d6663SPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess
4823c7d6663SPeter Brune - -snes_linesearch_max_it - The number of iterations for iterative line searches
483cd7522eaSPeter Brune 
484cd7522eaSPeter Brune    Notes:
485cd7522eaSPeter Brune    This is typically called from within a SNESSolve() implementation in order to
486cd7522eaSPeter Brune    help with convergence of the nonlinear method.  Various SNES types use line searches
487cd7522eaSPeter Brune    in different ways, but the overarching theme is that a line search is used to determine
488cd7522eaSPeter Brune    an optimal damping parameter of a step at each iteration of the method.  Each
489cd7522eaSPeter Brune    application of the line search may invoke SNESComputeFunction several times, and
490cd7522eaSPeter Brune    therefore may be fairly expensive.
491cd7522eaSPeter Brune 
492f40b411bSPeter Brune    Level: Intermediate
493f40b411bSPeter Brune 
494f1c6b773SPeter Brune .keywords: SNESLineSearch, Create
495f40b411bSPeter Brune 
496cd7522eaSPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchPreCheck(), SNESLineSearchPostCheck(), SNESSolve(), SNESComputeFunction()
497f40b411bSPeter Brune @*/
498bf388a1fSBarry Smith PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y)
499bf388a1fSBarry Smith {
500bf7f4e0aSPeter Brune   PetscErrorCode ierr;
501bf7f4e0aSPeter Brune 
502bf388a1fSBarry Smith   PetscFunctionBegin;
503f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
504bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
505bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
506bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
507bf7f4e0aSPeter Brune 
508bf7f4e0aSPeter Brune   linesearch->success = PETSC_TRUE;
509bf7f4e0aSPeter Brune 
510bf7f4e0aSPeter Brune   linesearch->vec_sol    = X;
511bf7f4e0aSPeter Brune   linesearch->vec_update = Y;
512bf7f4e0aSPeter Brune   linesearch->vec_func   = F;
513bf7f4e0aSPeter Brune 
514f1c6b773SPeter Brune   ierr = SNESLineSearchSetUp(linesearch);CHKERRQ(ierr);
515bf7f4e0aSPeter Brune 
516f5af7f23SKarl Rupp   if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */
517bf7f4e0aSPeter Brune 
518f5af7f23SKarl Rupp   if (fnorm) linesearch->fnorm = *fnorm;
519f5af7f23SKarl Rupp   else {
520bf7f4e0aSPeter Brune     ierr = VecNorm(F, NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
521bf7f4e0aSPeter Brune   }
522bf7f4e0aSPeter Brune 
523f1c6b773SPeter Brune   ierr = PetscLogEventBegin(SNESLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
524bf7f4e0aSPeter Brune 
525bf7f4e0aSPeter Brune   ierr = (*linesearch->ops->apply)(linesearch);CHKERRQ(ierr);
526bf7f4e0aSPeter Brune 
527f1c6b773SPeter Brune   ierr = PetscLogEventEnd(SNESLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
528bf7f4e0aSPeter Brune 
529f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
530bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
531bf7f4e0aSPeter Brune }
532bf7f4e0aSPeter Brune 
533bf7f4e0aSPeter Brune #undef __FUNCT__
534f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchDestroy"
535f40b411bSPeter Brune /*@
536f1c6b773SPeter Brune    SNESLineSearchDestroy - Destroys the line search instance.
537f40b411bSPeter Brune 
538f1c6b773SPeter Brune    Collective on SNESLineSearch
539f40b411bSPeter Brune 
540f40b411bSPeter Brune    Input Parameters:
5418e557f58SPeter Brune .  linesearch - The linesearch context
542f40b411bSPeter Brune 
543f40b411bSPeter Brune    Level: Intermediate
544f40b411bSPeter Brune 
54578bcb3b5SPeter Brune .keywords: SNESLineSearch, Destroy
546f40b411bSPeter Brune 
547cd7522eaSPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchReset(), SNESDestroy()
548f40b411bSPeter Brune @*/
549bf388a1fSBarry Smith PetscErrorCode SNESLineSearchDestroy(SNESLineSearch * linesearch)
550bf388a1fSBarry Smith {
551bf7f4e0aSPeter Brune   PetscErrorCode ierr;
552bf388a1fSBarry Smith 
553bf7f4e0aSPeter Brune   PetscFunctionBegin;
554bf7f4e0aSPeter Brune   if (!*linesearch) PetscFunctionReturn(0);
555f1c6b773SPeter Brune   PetscValidHeaderSpecific((*linesearch),SNESLINESEARCH_CLASSID,1);
556bf7f4e0aSPeter Brune   if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = 0; PetscFunctionReturn(0);}
557f05ece33SBarry Smith   ierr = PetscObjectAMSViewOff((PetscObject)*linesearch);CHKERRQ(ierr);
55822d28d08SBarry Smith   ierr = SNESLineSearchReset(*linesearch);CHKERRQ(ierr);
559f5af7f23SKarl Rupp   if ((*linesearch)->ops->destroy) (*linesearch)->ops->destroy(*linesearch);
560bf7f4e0aSPeter Brune   ierr = PetscViewerDestroy(&(*linesearch)->monitor);CHKERRQ(ierr);
561e7058c64SPeter Brune   ierr = PetscHeaderDestroy(linesearch);CHKERRQ(ierr);
562bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
563bf7f4e0aSPeter Brune }
564bf7f4e0aSPeter Brune 
565bf7f4e0aSPeter Brune #undef __FUNCT__
566f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetMonitor"
567f40b411bSPeter Brune /*@
568cd7522eaSPeter Brune    SNESLineSearchSetMonitor - Turns on/off printing useful information and debugging output about the line search.
569bf7f4e0aSPeter Brune 
570bf7f4e0aSPeter Brune    Input Parameters:
571bf7f4e0aSPeter Brune +  snes - nonlinear context obtained from SNESCreate()
572bf7f4e0aSPeter Brune -  flg - PETSC_TRUE to monitor the line search
573bf7f4e0aSPeter Brune 
574bf7f4e0aSPeter Brune    Logically Collective on SNES
575bf7f4e0aSPeter Brune 
576bf7f4e0aSPeter Brune    Options Database:
5778e557f58SPeter Brune .   -snes_linesearch_monitor - enables the monitor
578bf7f4e0aSPeter Brune 
579bf7f4e0aSPeter Brune    Level: intermediate
580bf7f4e0aSPeter Brune 
581bf7f4e0aSPeter Brune 
582cd7522eaSPeter Brune .seealso: SNESLineSearchGetMonitor(), PetscViewer
583bf7f4e0aSPeter Brune @*/
584f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetMonitor(SNESLineSearch linesearch, PetscBool flg)
585bf7f4e0aSPeter Brune {
586bf7f4e0aSPeter Brune   PetscErrorCode ierr;
587bf388a1fSBarry Smith 
588bf7f4e0aSPeter Brune   PetscFunctionBegin;
589bf7f4e0aSPeter Brune   if (flg && !linesearch->monitor) {
590ce94432eSBarry Smith     ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)linesearch),"stdout",&linesearch->monitor);CHKERRQ(ierr);
591bf7f4e0aSPeter Brune   } else if (!flg && linesearch->monitor) {
592bf7f4e0aSPeter Brune     ierr = PetscViewerDestroy(&linesearch->monitor);CHKERRQ(ierr);
593bf7f4e0aSPeter Brune   }
594bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
595bf7f4e0aSPeter Brune }
596bf7f4e0aSPeter Brune 
597bf7f4e0aSPeter Brune #undef __FUNCT__
598f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetMonitor"
599f40b411bSPeter Brune /*@
600cd7522eaSPeter Brune    SNESLineSearchGetMonitor - Gets the PetscViewer instance for the line search monitor.
6016a388c36SPeter Brune 
602f40b411bSPeter Brune    Input Parameters:
6038e557f58SPeter Brune .  linesearch - linesearch context
604f40b411bSPeter Brune 
605f40b411bSPeter Brune    Input Parameters:
6068e557f58SPeter Brune .  monitor - monitor context
607f40b411bSPeter Brune 
608f40b411bSPeter Brune    Logically Collective on SNES
609f40b411bSPeter Brune 
610f40b411bSPeter Brune 
611f40b411bSPeter Brune    Options Database Keys:
6128e557f58SPeter Brune .   -snes_linesearch_monitor - enables the monitor
613f40b411bSPeter Brune 
614f40b411bSPeter Brune    Level: intermediate
615f40b411bSPeter Brune 
616f40b411bSPeter Brune 
617cd7522eaSPeter Brune .seealso: SNESLineSearchSetMonitor(), PetscViewer
618f40b411bSPeter Brune @*/
619f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetMonitor(SNESLineSearch linesearch, PetscViewer *monitor)
6206a388c36SPeter Brune {
6216a388c36SPeter Brune   PetscFunctionBegin;
622f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
6236a388c36SPeter Brune   if (monitor) {
6246a388c36SPeter Brune     PetscValidPointer(monitor, 2);
6256a388c36SPeter Brune     *monitor = linesearch->monitor;
6266a388c36SPeter Brune   }
6276a388c36SPeter Brune   PetscFunctionReturn(0);
6286a388c36SPeter Brune }
6296a388c36SPeter Brune 
6306a388c36SPeter Brune #undef __FUNCT__
631f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetFromOptions"
632f40b411bSPeter Brune /*@
633f1c6b773SPeter Brune    SNESLineSearchSetFromOptions - Sets options for the line search
634f40b411bSPeter Brune 
635f40b411bSPeter Brune    Input Parameters:
6368e557f58SPeter Brune .  linesearch - linesearch context
637f40b411bSPeter Brune 
638f40b411bSPeter Brune    Options Database Keys:
6393c7d6663SPeter Brune + -snes_linesearch_type <type> - basic, bt, l2, cp, shell
6403c7d6663SPeter Brune . -snes_linesearch_order <order> - 1, 2, 3.  Most types only support certain orders (bt supports 2 or 3)
6415a9b6599SPeter Brune . -snes_linesearch_norms   - Turn on/off the linesearch norms for the basic linesearch type
64271eef1aeSPeter Brune . -snes_linesearch_minlambda - The minimum step length
6431a9b3a06SPeter Brune . -snes_linesearch_maxstep - The maximum step size
6441a9b3a06SPeter Brune . -snes_linesearch_rtol - Relative tolerance for iterative line searches
6451a9b3a06SPeter Brune . -snes_linesearch_atol - Absolute tolerance for iterative line searches
6461a9b3a06SPeter Brune . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches
6471a9b3a06SPeter Brune . -snes_linesearch_max_it - The number of iterations for iterative line searches
648cd7522eaSPeter Brune . -snes_linesearch_monitor - Print progress of line searches
6498e557f58SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter
650cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess.
6511a9b3a06SPeter Brune . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method
6521a9b3a06SPeter Brune - -snes_linesearch_precheck_picard_angle - Angle used in picard precheck method
653f40b411bSPeter Brune 
654f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
655f40b411bSPeter Brune 
656f40b411bSPeter Brune    Level: intermediate
657f40b411bSPeter Brune 
6583c7d6663SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetOrder(), SNESLineSearchSetType(), SNESLineSearchSetTolerances(), SNESLineSearchSetDamping(), SNESLineSearchPreCheckPicard()
659f40b411bSPeter Brune @*/
660bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch)
661bf388a1fSBarry Smith {
662bf7f4e0aSPeter Brune   PetscErrorCode ierr;
6631a4f838cSPeter Brune   const char     *deft = SNESLINESEARCHBASIC;
664bf7f4e0aSPeter Brune   char           type[256];
665bf7f4e0aSPeter Brune   PetscBool      flg, set;
666bf388a1fSBarry Smith 
667bf7f4e0aSPeter Brune   PetscFunctionBegin;
668607a6623SBarry Smith   if (!SNESLineSearchRegisterAllCalled) {ierr = SNESLineSearchRegisterAll();CHKERRQ(ierr);}
669bf7f4e0aSPeter Brune 
670bf7f4e0aSPeter Brune   ierr = PetscObjectOptionsBegin((PetscObject)linesearch);CHKERRQ(ierr);
671f5af7f23SKarl Rupp   if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name;
6723c7d6663SPeter Brune   ierr = PetscOptionsList("-snes_linesearch_type","Linesearch type","SNESLineSearchSetType",SNESLineSearchList,deft,type,256,&flg);CHKERRQ(ierr);
673bf7f4e0aSPeter Brune   if (flg) {
674f1c6b773SPeter Brune     ierr = SNESLineSearchSetType(linesearch,type);CHKERRQ(ierr);
675bf7f4e0aSPeter Brune   } else if (!((PetscObject)linesearch)->type_name) {
676f1c6b773SPeter Brune     ierr = SNESLineSearchSetType(linesearch,deft);CHKERRQ(ierr);
677bf7f4e0aSPeter Brune   }
678bf7f4e0aSPeter Brune 
6797a35526eSPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_monitor","Print progress of line searches","SNESSNESLineSearchSetMonitor",
680bf7f4e0aSPeter Brune                           linesearch->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
681f1c6b773SPeter Brune   if (set) {ierr = SNESLineSearchSetMonitor(linesearch,flg);CHKERRQ(ierr);}
682bf7f4e0aSPeter Brune 
6831a9b3a06SPeter Brune   /* tolerances */
68471eef1aeSPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_minlambda","Minimum step length","SNESLineSearchSetTolerances",linesearch->steptol,&linesearch->steptol,0);CHKERRQ(ierr);
6851a9b3a06SPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_maxstep","Maximum step size","SNESLineSearchSetTolerances",linesearch->maxstep,&linesearch->maxstep,0);CHKERRQ(ierr);
6861a9b3a06SPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_rtol","Relative tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->rtol,&linesearch->rtol,0);CHKERRQ(ierr);
6871a9b3a06SPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_atol","Absolute tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->atol,&linesearch->atol,0);CHKERRQ(ierr);
6881a9b3a06SPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_ltol","Change in lambda tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->ltol,&linesearch->ltol,0);CHKERRQ(ierr);
6898e557f58SPeter Brune   ierr = PetscOptionsInt("-snes_linesearch_max_it","Maximum iterations for iterative line searches","SNESLineSearchSetTolerances",linesearch->max_its,&linesearch->max_its,0);CHKERRQ(ierr);
6907a35526eSPeter Brune 
6911a9b3a06SPeter Brune   /* damping parameters */
6921a9b3a06SPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_damping","Line search damping and initial step guess","SNESLineSearchSetDamping",linesearch->damping,&linesearch->damping,0);CHKERRQ(ierr);
6931a9b3a06SPeter Brune 
6941a9b3a06SPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_keeplambda","Use previous lambda as damping","SNESLineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,0);CHKERRQ(ierr);
6951a9b3a06SPeter Brune 
6961a9b3a06SPeter Brune   /* precheck */
6977a35526eSPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_precheck_picard","Use a correction that sometimes improves convergence of Picard iteration","SNESLineSearchPreCheckPicard",flg,&flg,&set);CHKERRQ(ierr);
6987a35526eSPeter Brune   if (set) {
6997a35526eSPeter Brune     if (flg) {
7007a35526eSPeter Brune       linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */
701f5af7f23SKarl Rupp 
7027a35526eSPeter Brune       ierr = PetscOptionsReal("-snes_linesearch_precheck_picard_angle","Maximum angle at which to activate the correction",
7030298fd71SBarry Smith                               "none",linesearch->precheck_picard_angle,&linesearch->precheck_picard_angle,NULL);CHKERRQ(ierr);
7047a35526eSPeter Brune       ierr = SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&linesearch->precheck_picard_angle);CHKERRQ(ierr);
7057a35526eSPeter Brune     } else {
7060298fd71SBarry Smith       ierr = SNESLineSearchSetPreCheck(linesearch,NULL,NULL);CHKERRQ(ierr);
7077a35526eSPeter Brune     }
7087a35526eSPeter Brune   }
709b000cd8dSPeter Brune   ierr = PetscOptionsInt("-snes_linesearch_order","Order of approximation used in the line search","SNESLineSearchSetOrder",linesearch->order,&linesearch->order,0);CHKERRQ(ierr);
7101a9b3a06SPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_norms","Compute final norms in line search","SNESLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,0);CHKERRQ(ierr);
7117a35526eSPeter Brune 
7125a9b6599SPeter Brune   if (linesearch->ops->setfromoptions) {
7135a9b6599SPeter Brune     (*linesearch->ops->setfromoptions)(linesearch);CHKERRQ(ierr);
7145a9b6599SPeter Brune   }
7155a9b6599SPeter Brune 
716bf7f4e0aSPeter Brune   ierr = PetscObjectProcessOptionsHandlers((PetscObject)linesearch);CHKERRQ(ierr);
717bf7f4e0aSPeter Brune   ierr = PetscOptionsEnd();CHKERRQ(ierr);
718bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
719bf7f4e0aSPeter Brune }
720bf7f4e0aSPeter Brune 
721bf7f4e0aSPeter Brune #undef __FUNCT__
722f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchView"
723f40b411bSPeter Brune /*@
724cd7522eaSPeter Brune    SNESLineSearchView - Prints useful information about the line search not
725cd7522eaSPeter Brune    related to an individual call.
726f40b411bSPeter Brune 
727f40b411bSPeter Brune    Input Parameters:
7288e557f58SPeter Brune .  linesearch - linesearch context
729f40b411bSPeter Brune 
730f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
731f40b411bSPeter Brune 
732f40b411bSPeter Brune    Level: intermediate
733f40b411bSPeter Brune 
734f1c6b773SPeter Brune .seealso: SNESLineSearchCreate()
735f40b411bSPeter Brune @*/
736bf388a1fSBarry Smith PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer)
737bf388a1fSBarry Smith {
7387f1410a3SPeter Brune   PetscErrorCode ierr;
7397f1410a3SPeter Brune   PetscBool      iascii;
740bf388a1fSBarry Smith 
741bf7f4e0aSPeter Brune   PetscFunctionBegin;
7427f1410a3SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
7437f1410a3SPeter Brune   if (!viewer) {
744ce94432eSBarry Smith     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch),&viewer);CHKERRQ(ierr);
7457f1410a3SPeter Brune   }
7467f1410a3SPeter Brune   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
7477f1410a3SPeter Brune   PetscCheckSameComm(linesearch,1,viewer,2);
748f40b411bSPeter Brune 
749251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
7507f1410a3SPeter Brune   if (iascii) {
7517f1410a3SPeter Brune     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)linesearch,viewer,"SNESLineSearch Object");CHKERRQ(ierr);
7527f1410a3SPeter Brune     if (linesearch->ops->view) {
7537f1410a3SPeter Brune       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
7547f1410a3SPeter Brune       ierr = (*linesearch->ops->view)(linesearch,viewer);CHKERRQ(ierr);
7557f1410a3SPeter Brune       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
7567f1410a3SPeter Brune     }
757c69d1a72SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep,(double)linesearch->steptol);CHKERRQ(ierr);
758c69d1a72SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol,(double)linesearch->atol,(double)linesearch->ltol);CHKERRQ(ierr);
7597f1410a3SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D\n", linesearch->max_its);CHKERRQ(ierr);
7606b2b7091SBarry Smith     if (linesearch->ops->precheck) {
7616b2b7091SBarry Smith       if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) {
7627f1410a3SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  using precheck step to speed up Picard convergence\n", linesearch->max_its);CHKERRQ(ierr);
7637f1410a3SPeter Brune       } else {
7647f1410a3SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  using user-defined precheck step\n", linesearch->max_its);CHKERRQ(ierr);
7657f1410a3SPeter Brune       }
7667f1410a3SPeter Brune     }
7676b2b7091SBarry Smith     if (linesearch->ops->postcheck) {
7687f1410a3SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer,"  using user-defined postcheck step\n", linesearch->max_its);CHKERRQ(ierr);
7697f1410a3SPeter Brune     }
7707f1410a3SPeter Brune   }
771bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
772bf7f4e0aSPeter Brune }
773bf7f4e0aSPeter Brune 
774bf7f4e0aSPeter Brune #undef __FUNCT__
775f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetType"
776ea5d4fccSPeter Brune /*@C
777f1c6b773SPeter Brune    SNESLineSearchSetType - Sets the linesearch type
778f40b411bSPeter Brune 
779f40b411bSPeter Brune    Input Parameters:
7808e557f58SPeter Brune +  linesearch - linesearch context
781f40b411bSPeter Brune -  type - The type of line search to be used
782f40b411bSPeter Brune 
783cd7522eaSPeter Brune    Available Types:
784cd7522eaSPeter Brune +  basic - Simple damping line search.
7858e557f58SPeter Brune .  bt - Backtracking line search over the L2 norm of the function
7868e557f58SPeter Brune .  l2 - Secant line search over the L2 norm of the function
7878e557f58SPeter Brune .  cp - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x)
7888e557f58SPeter Brune -  shell - User provided SNESLineSearch implementation
789cd7522eaSPeter Brune 
790f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
791f40b411bSPeter Brune 
792f40b411bSPeter Brune    Level: intermediate
793f40b411bSPeter Brune 
794f40b411bSPeter Brune 
795f1c6b773SPeter Brune .seealso: SNESLineSearchCreate()
796f40b411bSPeter Brune @*/
79719fd82e9SBarry Smith PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type)
798bf7f4e0aSPeter Brune {
799f1c6b773SPeter Brune   PetscErrorCode ierr,(*r)(SNESLineSearch);
800bf7f4e0aSPeter Brune   PetscBool      match;
801bf7f4e0aSPeter Brune 
802bf7f4e0aSPeter Brune   PetscFunctionBegin;
803f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
804bf7f4e0aSPeter Brune   PetscValidCharPointer(type,2);
805bf7f4e0aSPeter Brune 
806251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)linesearch,type,&match);CHKERRQ(ierr);
807bf7f4e0aSPeter Brune   if (match) PetscFunctionReturn(0);
808bf7f4e0aSPeter Brune 
8091c9cd337SJed Brown   ierr = PetscFunctionListFind(SNESLineSearchList,type,&r);CHKERRQ(ierr);
810bf7f4e0aSPeter Brune   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type);
811bf7f4e0aSPeter Brune   /* Destroy the previous private linesearch context */
812bf7f4e0aSPeter Brune   if (linesearch->ops->destroy) {
813bf7f4e0aSPeter Brune     ierr = (*(linesearch)->ops->destroy)(linesearch);CHKERRQ(ierr);
814f5af7f23SKarl Rupp 
8150298fd71SBarry Smith     linesearch->ops->destroy = NULL;
816bf7f4e0aSPeter Brune   }
817f1c6b773SPeter Brune   /* Reinitialize function pointers in SNESLineSearchOps structure */
818bf7f4e0aSPeter Brune   linesearch->ops->apply          = 0;
819bf7f4e0aSPeter Brune   linesearch->ops->view           = 0;
820bf7f4e0aSPeter Brune   linesearch->ops->setfromoptions = 0;
821bf7f4e0aSPeter Brune   linesearch->ops->destroy        = 0;
822bf7f4e0aSPeter Brune 
823bf7f4e0aSPeter Brune   ierr = PetscObjectChangeTypeName((PetscObject)linesearch,type);CHKERRQ(ierr);
824bf7f4e0aSPeter Brune   ierr = (*r)(linesearch);CHKERRQ(ierr);
825bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
826bf7f4e0aSPeter Brune }
827bf7f4e0aSPeter Brune 
828bf7f4e0aSPeter Brune #undef __FUNCT__
829f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetSNES"
830f40b411bSPeter Brune /*@
83178bcb3b5SPeter Brune    SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation.
832f40b411bSPeter Brune 
833f40b411bSPeter Brune    Input Parameters:
8348e557f58SPeter Brune +  linesearch - linesearch context
835f40b411bSPeter Brune -  snes - The snes instance
836f40b411bSPeter Brune 
83778bcb3b5SPeter Brune    Level: developer
83878bcb3b5SPeter Brune 
83978bcb3b5SPeter Brune    Notes:
84078bcb3b5SPeter Brune    This happens automatically when the line search is gotten/created with
841*7601faf0SJed Brown    SNESGetLineSearch().  This routine is therefore mainly called within SNES
84278bcb3b5SPeter Brune    implementations.
843f40b411bSPeter Brune 
8448141a3b9SPeter Brune    Level: developer
845f40b411bSPeter Brune 
846cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
847f40b411bSPeter Brune @*/
848bf388a1fSBarry Smith PetscErrorCode  SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes)
849bf388a1fSBarry Smith {
850bf7f4e0aSPeter Brune   PetscFunctionBegin;
851f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
852bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,2);
853bf7f4e0aSPeter Brune   linesearch->snes = snes;
854bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
855bf7f4e0aSPeter Brune }
856bf7f4e0aSPeter Brune 
857bf7f4e0aSPeter Brune #undef __FUNCT__
858f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetSNES"
859f40b411bSPeter Brune /*@
8608141a3b9SPeter Brune    SNESLineSearchGetSNES - Gets the SNES instance associated with the line search.
8618141a3b9SPeter Brune    Having an associated SNES is necessary because most line search implementations must be able to
8628141a3b9SPeter Brune    evaluate the function using SNESComputeFunction() for the associated SNES.  This routine
8638141a3b9SPeter Brune    is used in the line search implementations when one must get this associated SNES instance.
864f40b411bSPeter Brune 
865f40b411bSPeter Brune    Input Parameters:
8668e557f58SPeter Brune .  linesearch - linesearch context
867f40b411bSPeter Brune 
868f40b411bSPeter Brune    Output Parameters:
869f40b411bSPeter Brune .  snes - The snes instance
870f40b411bSPeter Brune 
8718141a3b9SPeter Brune    Level: developer
872f40b411bSPeter Brune 
873cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
874f40b411bSPeter Brune @*/
875bf388a1fSBarry Smith PetscErrorCode  SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes)
876bf388a1fSBarry Smith {
877bf7f4e0aSPeter Brune   PetscFunctionBegin;
878f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
8796a388c36SPeter Brune   PetscValidPointer(snes, 2);
880bf7f4e0aSPeter Brune   *snes = linesearch->snes;
881bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
882bf7f4e0aSPeter Brune }
883bf7f4e0aSPeter Brune 
8846a388c36SPeter Brune #undef __FUNCT__
885f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetLambda"
886f40b411bSPeter Brune /*@
887f1c6b773SPeter Brune    SNESLineSearchGetLambda - Gets the last linesearch steplength discovered.
888f40b411bSPeter Brune 
889f40b411bSPeter Brune    Input Parameters:
8908e557f58SPeter Brune .  linesearch - linesearch context
891f40b411bSPeter Brune 
892f40b411bSPeter Brune    Output Parameters:
893cd7522eaSPeter Brune .  lambda - The last steplength computed during SNESLineSearchApply()
894f40b411bSPeter Brune 
89578bcb3b5SPeter Brune    Level: advanced
89678bcb3b5SPeter Brune 
8978e557f58SPeter Brune    Notes:
8988e557f58SPeter Brune    This is useful in methods where the solver is ill-scaled and
89978bcb3b5SPeter Brune    requires some adaptive notion of the difference in scale between the
90078bcb3b5SPeter Brune    solution and the function.  For instance, SNESQN may be scaled by the
90178bcb3b5SPeter Brune    line search lambda using the argument -snes_qn_scaling ls.
90278bcb3b5SPeter Brune 
903f40b411bSPeter Brune 
904cd7522eaSPeter Brune .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetDamping(), SNESLineSearchApply()
905f40b411bSPeter Brune @*/
906f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda)
9076a388c36SPeter Brune {
9086a388c36SPeter Brune   PetscFunctionBegin;
909f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
9106a388c36SPeter Brune   PetscValidPointer(lambda, 2);
9116a388c36SPeter Brune   *lambda = linesearch->lambda;
9126a388c36SPeter Brune   PetscFunctionReturn(0);
9136a388c36SPeter Brune }
9146a388c36SPeter Brune 
9156a388c36SPeter Brune #undef __FUNCT__
916f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetLambda"
917f40b411bSPeter Brune /*@
918f1c6b773SPeter Brune    SNESLineSearchSetLambda - Sets the linesearch steplength.
919f40b411bSPeter Brune 
920f40b411bSPeter Brune    Input Parameters:
9218e557f58SPeter Brune +  linesearch - linesearch context
922f40b411bSPeter Brune -  lambda - The last steplength.
923f40b411bSPeter Brune 
924cd7522eaSPeter Brune    Notes:
925cd7522eaSPeter Brune    This routine is typically used within implementations of SNESLineSearchApply
926cd7522eaSPeter Brune    to set the final steplength.  This routine (and SNESLineSearchGetLambda()) were
927cd7522eaSPeter Brune    added in order to facilitate Quasi-Newton methods that use the previous steplength
928cd7522eaSPeter Brune    as an inner scaling parameter.
929cd7522eaSPeter Brune 
93078bcb3b5SPeter Brune    Level: advanced
931f40b411bSPeter Brune 
932f1c6b773SPeter Brune .seealso: SNESLineSearchGetLambda()
933f40b411bSPeter Brune @*/
934f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
9356a388c36SPeter Brune {
9366a388c36SPeter Brune   PetscFunctionBegin;
937f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
9386a388c36SPeter Brune   linesearch->lambda = lambda;
9396a388c36SPeter Brune   PetscFunctionReturn(0);
9406a388c36SPeter Brune }
9416a388c36SPeter Brune 
9426a388c36SPeter Brune #undef  __FUNCT__
943f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetTolerances"
944f40b411bSPeter Brune /*@
9453c7d6663SPeter Brune    SNESLineSearchGetTolerances - Gets the tolerances for the linesearch.  These include
94678bcb3b5SPeter Brune    tolerances for the relative and absolute change in the function norm, the change
94778bcb3b5SPeter Brune    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
94878bcb3b5SPeter Brune    and the maximum number of iterations the line search procedure may take.
949f40b411bSPeter Brune 
950f40b411bSPeter Brune    Input Parameters:
9518e557f58SPeter Brune .  linesearch - linesearch context
952f40b411bSPeter Brune 
953f40b411bSPeter Brune    Output Parameters:
954516fe3c3SPeter Brune +  steptol - The minimum steplength
9556cc8e53bSPeter Brune .  maxstep - The maximum steplength
956516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
957516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
958516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
959516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
960f40b411bSPeter Brune 
96178bcb3b5SPeter Brune    Level: intermediate
96278bcb3b5SPeter Brune 
96378bcb3b5SPeter Brune    Notes:
96478bcb3b5SPeter Brune    Different line searches may implement these parameters slightly differently as
9653c7d6663SPeter Brune    the type requires.
966516fe3c3SPeter Brune 
967f1c6b773SPeter Brune .seealso: SNESLineSearchSetTolerances()
968f40b411bSPeter Brune @*/
969f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
9706a388c36SPeter Brune {
9716a388c36SPeter Brune   PetscFunctionBegin;
972f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
973516fe3c3SPeter Brune   if (steptol) {
9746a388c36SPeter Brune     PetscValidPointer(steptol, 2);
9756a388c36SPeter Brune     *steptol = linesearch->steptol;
976516fe3c3SPeter Brune   }
977516fe3c3SPeter Brune   if (maxstep) {
978516fe3c3SPeter Brune     PetscValidPointer(maxstep, 3);
979516fe3c3SPeter Brune     *maxstep = linesearch->maxstep;
980516fe3c3SPeter Brune   }
981516fe3c3SPeter Brune   if (rtol) {
982516fe3c3SPeter Brune     PetscValidPointer(rtol, 4);
983516fe3c3SPeter Brune     *rtol = linesearch->rtol;
984516fe3c3SPeter Brune   }
985516fe3c3SPeter Brune   if (atol) {
986516fe3c3SPeter Brune     PetscValidPointer(atol, 5);
987516fe3c3SPeter Brune     *atol = linesearch->atol;
988516fe3c3SPeter Brune   }
989516fe3c3SPeter Brune   if (ltol) {
990516fe3c3SPeter Brune     PetscValidPointer(ltol, 6);
991516fe3c3SPeter Brune     *ltol = linesearch->ltol;
992516fe3c3SPeter Brune   }
993516fe3c3SPeter Brune   if (max_its) {
994516fe3c3SPeter Brune     PetscValidPointer(max_its, 7);
995516fe3c3SPeter Brune     *max_its = linesearch->max_its;
996516fe3c3SPeter Brune   }
9976a388c36SPeter Brune   PetscFunctionReturn(0);
9986a388c36SPeter Brune }
9996a388c36SPeter Brune 
10006a388c36SPeter Brune #undef  __FUNCT__
1001f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetTolerances"
1002f40b411bSPeter Brune /*@
10033c7d6663SPeter Brune    SNESLineSearchSetTolerances -  Gets the tolerances for the linesearch.  These include
100478bcb3b5SPeter Brune    tolerances for the relative and absolute change in the function norm, the change
100578bcb3b5SPeter Brune    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
100678bcb3b5SPeter Brune    and the maximum number of iterations the line search procedure may take.
1007f40b411bSPeter Brune 
1008f40b411bSPeter Brune    Input Parameters:
10098e557f58SPeter Brune +  linesearch - linesearch context
1010516fe3c3SPeter Brune .  steptol - The minimum steplength
10116cc8e53bSPeter Brune .  maxstep - The maximum steplength
1012516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
1013516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
1014516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
1015516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
1016f40b411bSPeter Brune 
101778bcb3b5SPeter Brune    Notes:
10183c7d6663SPeter Brune    The user may choose to not set any of the tolerances using PETSC_DEFAULT in
101978bcb3b5SPeter Brune    place of an argument.
1020f40b411bSPeter Brune 
102178bcb3b5SPeter Brune    Level: intermediate
1022516fe3c3SPeter Brune 
1023f1c6b773SPeter Brune .seealso: SNESLineSearchGetTolerances()
1024f40b411bSPeter Brune @*/
1025f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its)
10266a388c36SPeter Brune {
10276a388c36SPeter Brune   PetscFunctionBegin;
1028f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1029d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,steptol,2);
1030d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,maxstep,3);
1031d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,rtol,4);
1032d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,atol,5);
1033d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,ltol,6);
1034d3952184SSatish Balay   PetscValidLogicalCollectiveInt(linesearch,max_its,7);
1035d3952184SSatish Balay 
1036d3952184SSatish Balay   if (steptol!= PETSC_DEFAULT) {
1037ce94432eSBarry Smith     if (steptol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Minimum step length %14.12e must be non-negative",(double)steptol);
10386a388c36SPeter Brune     linesearch->steptol = steptol;
1039d3952184SSatish Balay   }
1040d3952184SSatish Balay 
1041d3952184SSatish Balay   if (maxstep!= PETSC_DEFAULT) {
1042ce94432eSBarry Smith     if (maxstep < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum step length %14.12e must be non-negative",(double)maxstep);
1043516fe3c3SPeter Brune     linesearch->maxstep = maxstep;
1044d3952184SSatish Balay   }
1045d3952184SSatish Balay 
1046d3952184SSatish Balay   if (rtol != PETSC_DEFAULT) {
1047ce94432eSBarry 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);
1048516fe3c3SPeter Brune     linesearch->rtol = rtol;
1049d3952184SSatish Balay   }
1050d3952184SSatish Balay 
1051d3952184SSatish Balay   if (atol != PETSC_DEFAULT) {
1052ce94432eSBarry Smith     if (atol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %14.12e must be non-negative",(double)atol);
1053516fe3c3SPeter Brune     linesearch->atol = atol;
1054d3952184SSatish Balay   }
1055d3952184SSatish Balay 
1056d3952184SSatish Balay   if (ltol != PETSC_DEFAULT) {
1057ce94432eSBarry Smith     if (ltol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Labmda tolerance %14.12e must be non-negative",(double)ltol);
1058516fe3c3SPeter Brune     linesearch->ltol = ltol;
1059d3952184SSatish Balay   }
1060d3952184SSatish Balay 
1061d3952184SSatish Balay   if (max_its != PETSC_DEFAULT) {
1062ce94432eSBarry Smith     if (max_its < 0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",max_its);
1063516fe3c3SPeter Brune     linesearch->max_its = max_its;
1064d3952184SSatish Balay   }
10656a388c36SPeter Brune   PetscFunctionReturn(0);
10666a388c36SPeter Brune }
10676a388c36SPeter Brune 
10686a388c36SPeter Brune #undef __FUNCT__
1069f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetDamping"
1070f40b411bSPeter Brune /*@
1071f1c6b773SPeter Brune    SNESLineSearchGetDamping - Gets the line search damping parameter.
1072f40b411bSPeter Brune 
1073f40b411bSPeter Brune    Input Parameters:
10748e557f58SPeter Brune .  linesearch - linesearch context
1075f40b411bSPeter Brune 
1076f40b411bSPeter Brune    Output Parameters:
10778e557f58SPeter Brune .  damping - The damping parameter
1078f40b411bSPeter Brune 
107978bcb3b5SPeter Brune    Level: advanced
1080f40b411bSPeter Brune 
108178bcb3b5SPeter Brune .seealso: SNESLineSearchGetStepTolerance(), SNESQN
1082f40b411bSPeter Brune @*/
1083f40b411bSPeter Brune 
1084f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping)
10856a388c36SPeter Brune {
10866a388c36SPeter Brune   PetscFunctionBegin;
1087f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
10886a388c36SPeter Brune   PetscValidPointer(damping, 2);
10896a388c36SPeter Brune   *damping = linesearch->damping;
10906a388c36SPeter Brune   PetscFunctionReturn(0);
10916a388c36SPeter Brune }
10926a388c36SPeter Brune 
10936a388c36SPeter Brune #undef __FUNCT__
1094f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetDamping"
1095f40b411bSPeter Brune /*@
1096f1c6b773SPeter Brune    SNESLineSearchSetDamping - Sets the line search damping paramter.
1097f40b411bSPeter Brune 
1098f40b411bSPeter Brune    Input Parameters:
109978bcb3b5SPeter Brune .  linesearch - linesearch context
110078bcb3b5SPeter Brune .  damping - The damping parameter
1101f40b411bSPeter Brune 
1102f40b411bSPeter Brune    Level: intermediate
1103f40b411bSPeter Brune 
1104cd7522eaSPeter Brune    Notes:
1105cd7522eaSPeter Brune    The basic line search merely takes the update step scaled by the damping parameter.
1106cd7522eaSPeter Brune    The use of the damping parameter in the l2 and cp line searches is much more subtle;
110778bcb3b5SPeter Brune    it is used as a starting point in calculating the secant step. However, the eventual
1108cd7522eaSPeter Brune    step may be of greater length than the damping parameter.  In the bt line search it is
1109cd7522eaSPeter Brune    used as the maximum possible step length, as the bt line search only backtracks.
1110cd7522eaSPeter Brune 
1111f1c6b773SPeter Brune .seealso: SNESLineSearchGetDamping()
1112f40b411bSPeter Brune @*/
1113f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping)
11146a388c36SPeter Brune {
11156a388c36SPeter Brune   PetscFunctionBegin;
1116f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
11176a388c36SPeter Brune   linesearch->damping = damping;
11186a388c36SPeter Brune   PetscFunctionReturn(0);
11196a388c36SPeter Brune }
11206a388c36SPeter Brune 
11216a388c36SPeter Brune #undef __FUNCT__
112259405d5eSPeter Brune #define __FUNCT__ "SNESLineSearchGetOrder"
112359405d5eSPeter Brune /*@
112459405d5eSPeter Brune    SNESLineSearchGetOrder - Gets the line search approximation order.
112559405d5eSPeter Brune 
112659405d5eSPeter Brune    Input Parameters:
112778bcb3b5SPeter Brune .  linesearch - linesearch context
112859405d5eSPeter Brune 
112959405d5eSPeter Brune    Output Parameters:
11308e557f58SPeter Brune .  order - The order
113159405d5eSPeter Brune 
113278bcb3b5SPeter Brune    Possible Values for order:
11333c7d6663SPeter Brune +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
11343c7d6663SPeter Brune .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
11353c7d6663SPeter Brune -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
113678bcb3b5SPeter Brune 
113759405d5eSPeter Brune    Level: intermediate
113859405d5eSPeter Brune 
113959405d5eSPeter Brune .seealso: SNESLineSearchSetOrder()
114059405d5eSPeter Brune @*/
114159405d5eSPeter Brune 
1142b000cd8dSPeter Brune PetscErrorCode  SNESLineSearchGetOrder(SNESLineSearch linesearch,PetscInt *order)
114359405d5eSPeter Brune {
114459405d5eSPeter Brune   PetscFunctionBegin;
114559405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
114659405d5eSPeter Brune   PetscValidPointer(order, 2);
114759405d5eSPeter Brune   *order = linesearch->order;
114859405d5eSPeter Brune   PetscFunctionReturn(0);
114959405d5eSPeter Brune }
115059405d5eSPeter Brune 
115159405d5eSPeter Brune #undef __FUNCT__
115259405d5eSPeter Brune #define __FUNCT__ "SNESLineSearchSetOrder"
115359405d5eSPeter Brune /*@
115459405d5eSPeter Brune    SNESLineSearchSetOrder - Sets the line search damping paramter.
115559405d5eSPeter Brune 
115659405d5eSPeter Brune    Input Parameters:
115778bcb3b5SPeter Brune .  linesearch - linesearch context
115878bcb3b5SPeter Brune .  order - The damping parameter
115959405d5eSPeter Brune 
116059405d5eSPeter Brune    Level: intermediate
116159405d5eSPeter Brune 
116278bcb3b5SPeter Brune    Possible Values for order:
11633c7d6663SPeter Brune +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
11643c7d6663SPeter Brune .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
11653c7d6663SPeter Brune -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
116678bcb3b5SPeter Brune 
116759405d5eSPeter Brune    Notes:
116859405d5eSPeter Brune    Variable orders are supported by the following line searches:
116978bcb3b5SPeter Brune +  bt - cubic and quadratic
117078bcb3b5SPeter Brune -  cp - linear and quadratic
117159405d5eSPeter Brune 
117259405d5eSPeter Brune .seealso: SNESLineSearchGetOrder()
117359405d5eSPeter Brune @*/
1174b000cd8dSPeter Brune PetscErrorCode  SNESLineSearchSetOrder(SNESLineSearch linesearch,PetscInt order)
117559405d5eSPeter Brune {
117659405d5eSPeter Brune   PetscFunctionBegin;
117759405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
117859405d5eSPeter Brune   linesearch->order = order;
117959405d5eSPeter Brune   PetscFunctionReturn(0);
118059405d5eSPeter Brune }
118159405d5eSPeter Brune 
118259405d5eSPeter Brune #undef __FUNCT__
1183f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetNorms"
1184f40b411bSPeter Brune /*@
1185f1c6b773SPeter Brune    SNESLineSearchGetNorms - Gets the norms for for X, Y, and F.
1186f40b411bSPeter Brune 
1187f40b411bSPeter Brune    Input Parameters:
118878bcb3b5SPeter Brune .  linesearch - linesearch context
1189f40b411bSPeter Brune 
1190f40b411bSPeter Brune    Output Parameters:
1191f40b411bSPeter Brune +  xnorm - The norm of the current solution
1192f40b411bSPeter Brune .  fnorm - The norm of the current function
1193f40b411bSPeter Brune -  ynorm - The norm of the current update
1194f40b411bSPeter Brune 
1195cd7522eaSPeter Brune    Notes:
1196cd7522eaSPeter Brune    This function is mainly called from SNES implementations.
1197cd7522eaSPeter Brune 
119878bcb3b5SPeter Brune    Level: developer
1199f40b411bSPeter Brune 
1200f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms() SNESLineSearchGetVecs()
1201f40b411bSPeter Brune @*/
1202f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm)
1203bf7f4e0aSPeter Brune {
1204bf7f4e0aSPeter Brune   PetscFunctionBegin;
1205f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1206f5af7f23SKarl Rupp   if (xnorm) *xnorm = linesearch->xnorm;
1207f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
1208f5af7f23SKarl Rupp   if (ynorm) *ynorm = linesearch->ynorm;
1209bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1210bf7f4e0aSPeter Brune }
1211bf7f4e0aSPeter Brune 
1212e7058c64SPeter Brune #undef __FUNCT__
1213f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetNorms"
1214f40b411bSPeter Brune /*@
1215f1c6b773SPeter Brune    SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F.
1216f40b411bSPeter Brune 
1217f40b411bSPeter Brune    Input Parameters:
121878bcb3b5SPeter Brune +  linesearch - linesearch context
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 
122378bcb3b5SPeter Brune    Level: advanced
1224f40b411bSPeter Brune 
1225f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1226f40b411bSPeter Brune @*/
1227f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
12286a388c36SPeter Brune {
12296a388c36SPeter Brune   PetscFunctionBegin;
1230f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
12316a388c36SPeter Brune   linesearch->xnorm = xnorm;
12326a388c36SPeter Brune   linesearch->fnorm = fnorm;
12336a388c36SPeter Brune   linesearch->ynorm = ynorm;
12346a388c36SPeter Brune   PetscFunctionReturn(0);
12356a388c36SPeter Brune }
12366a388c36SPeter Brune 
12376a388c36SPeter Brune #undef __FUNCT__
1238f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchComputeNorms"
1239f40b411bSPeter Brune /*@
1240f1c6b773SPeter Brune    SNESLineSearchComputeNorms - Computes the norms of X, F, and Y.
1241f40b411bSPeter Brune 
1242f40b411bSPeter Brune    Input Parameters:
124378bcb3b5SPeter Brune .  linesearch - linesearch context
1244f40b411bSPeter Brune 
1245f40b411bSPeter Brune    Options Database Keys:
12468e557f58SPeter Brune .   -snes_linesearch_norms - turn norm computation on or off
1247f40b411bSPeter Brune 
1248f40b411bSPeter Brune    Level: intermediate
1249f40b411bSPeter Brune 
125078bcb3b5SPeter Brune .seealso: SNESLineSearchGetNorms, SNESLineSearchSetNorms(), SNESLineSearchSetComputeNorms()
1251f40b411bSPeter Brune @*/
1252f1c6b773SPeter Brune PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
12536a388c36SPeter Brune {
12546a388c36SPeter Brune   PetscErrorCode ierr;
12559bd66eb0SPeter Brune   SNES           snes;
1256bf388a1fSBarry Smith 
12576a388c36SPeter Brune   PetscFunctionBegin;
12586a388c36SPeter Brune   if (linesearch->norms) {
12599bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
1260f1c6b773SPeter Brune       ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
12619bd66eb0SPeter Brune       ierr = VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
12629bd66eb0SPeter Brune       ierr = VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
12639bd66eb0SPeter Brune       ierr = (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);CHKERRQ(ierr);
12649bd66eb0SPeter Brune     } else {
12656a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_func,   NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
12666a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_sol,    NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
12676a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
12686a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_func,     NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
12696a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_sol,      NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
12706a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_update,   NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
12716a388c36SPeter Brune     }
12729bd66eb0SPeter Brune   }
12736a388c36SPeter Brune   PetscFunctionReturn(0);
12746a388c36SPeter Brune }
12756a388c36SPeter Brune 
12766f263ca3SPeter Brune #undef __FUNCT__
12776f263ca3SPeter Brune #define __FUNCT__ "SNESLineSearchSetComputeNorms"
12786f263ca3SPeter Brune /*@
12796f263ca3SPeter Brune    SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search.
12806f263ca3SPeter Brune 
12816f263ca3SPeter Brune    Input Parameters:
128278bcb3b5SPeter Brune +  linesearch  - linesearch context
128378bcb3b5SPeter Brune -  flg  - indicates whether or not to compute norms
12846f263ca3SPeter Brune 
12856f263ca3SPeter Brune    Options Database Keys:
12868e557f58SPeter Brune .   -snes_linesearch_norms - turn norm computation on or off
12876f263ca3SPeter Brune 
12886f263ca3SPeter Brune    Notes:
12891a4f838cSPeter Brune    This is most relevant to the SNESLINESEARCHBASIC line search type.
12906f263ca3SPeter Brune 
12916f263ca3SPeter Brune    Level: intermediate
12926f263ca3SPeter Brune 
12931a4f838cSPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetNorms(), SNESLineSearchComputeNorms(), SNESLINESEARCHBASIC
12946f263ca3SPeter Brune @*/
12956f263ca3SPeter Brune PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg)
12966f263ca3SPeter Brune {
12976f263ca3SPeter Brune   PetscFunctionBegin;
12986f263ca3SPeter Brune   linesearch->norms = flg;
12996f263ca3SPeter Brune   PetscFunctionReturn(0);
13006f263ca3SPeter Brune }
13016f263ca3SPeter Brune 
13026a388c36SPeter Brune #undef __FUNCT__
1303f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetVecs"
1304f40b411bSPeter Brune /*@
1305f1c6b773SPeter Brune    SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context
1306f40b411bSPeter Brune 
1307f40b411bSPeter Brune    Input Parameters:
130878bcb3b5SPeter Brune .  linesearch - linesearch context
1309f40b411bSPeter Brune 
1310f40b411bSPeter Brune    Output Parameters:
1311f40b411bSPeter Brune +  X - The old solution
1312f40b411bSPeter Brune .  F - The old function
1313f40b411bSPeter Brune .  Y - The search direction
1314f40b411bSPeter Brune .  W - The new solution
1315f40b411bSPeter Brune -  G - The new function
1316f40b411bSPeter Brune 
131778bcb3b5SPeter Brune    Level: advanced
1318f40b411bSPeter Brune 
1319f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1320f40b411bSPeter Brune @*/
1321bf388a1fSBarry Smith PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G)
1322bf388a1fSBarry Smith {
13236a388c36SPeter Brune   PetscFunctionBegin;
1324f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
13256a388c36SPeter Brune   if (X) {
13266a388c36SPeter Brune     PetscValidPointer(X, 2);
13276a388c36SPeter Brune     *X = linesearch->vec_sol;
13286a388c36SPeter Brune   }
13296a388c36SPeter Brune   if (F) {
13306a388c36SPeter Brune     PetscValidPointer(F, 3);
13316a388c36SPeter Brune     *F = linesearch->vec_func;
13326a388c36SPeter Brune   }
13336a388c36SPeter Brune   if (Y) {
13346a388c36SPeter Brune     PetscValidPointer(Y, 4);
13356a388c36SPeter Brune     *Y = linesearch->vec_update;
13366a388c36SPeter Brune   }
13376a388c36SPeter Brune   if (W) {
13386a388c36SPeter Brune     PetscValidPointer(W, 5);
13396a388c36SPeter Brune     *W = linesearch->vec_sol_new;
13406a388c36SPeter Brune   }
13416a388c36SPeter Brune   if (G) {
13426a388c36SPeter Brune     PetscValidPointer(G, 6);
13436a388c36SPeter Brune     *G = linesearch->vec_func_new;
13446a388c36SPeter Brune   }
13456a388c36SPeter Brune   PetscFunctionReturn(0);
13466a388c36SPeter Brune }
13476a388c36SPeter Brune 
13486a388c36SPeter Brune #undef __FUNCT__
1349f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetVecs"
1350f40b411bSPeter Brune /*@
1351f1c6b773SPeter Brune    SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context
1352f40b411bSPeter Brune 
1353f40b411bSPeter Brune    Input Parameters:
135478bcb3b5SPeter Brune +  linesearch - linesearch context
1355f40b411bSPeter Brune .  X - The old solution
1356f40b411bSPeter Brune .  F - The old function
1357f40b411bSPeter Brune .  Y - The search direction
1358f40b411bSPeter Brune .  W - The new solution
1359f40b411bSPeter Brune -  G - The new function
1360f40b411bSPeter Brune 
136178bcb3b5SPeter Brune    Level: advanced
1362f40b411bSPeter Brune 
1363f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms(), SNESLineSearchGetVecs()
1364f40b411bSPeter Brune @*/
1365bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G)
1366bf388a1fSBarry Smith {
13676a388c36SPeter Brune   PetscFunctionBegin;
1368f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
13696a388c36SPeter Brune   if (X) {
13706a388c36SPeter Brune     PetscValidHeaderSpecific(X,VEC_CLASSID,2);
13716a388c36SPeter Brune     linesearch->vec_sol = X;
13726a388c36SPeter Brune   }
13736a388c36SPeter Brune   if (F) {
13746a388c36SPeter Brune     PetscValidHeaderSpecific(F,VEC_CLASSID,3);
13756a388c36SPeter Brune     linesearch->vec_func = F;
13766a388c36SPeter Brune   }
13776a388c36SPeter Brune   if (Y) {
13786a388c36SPeter Brune     PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
13796a388c36SPeter Brune     linesearch->vec_update = Y;
13806a388c36SPeter Brune   }
13816a388c36SPeter Brune   if (W) {
13826a388c36SPeter Brune     PetscValidHeaderSpecific(W,VEC_CLASSID,5);
13836a388c36SPeter Brune     linesearch->vec_sol_new = W;
13846a388c36SPeter Brune   }
13856a388c36SPeter Brune   if (G) {
13866a388c36SPeter Brune     PetscValidHeaderSpecific(G,VEC_CLASSID,6);
13876a388c36SPeter Brune     linesearch->vec_func_new = G;
13886a388c36SPeter Brune   }
13896a388c36SPeter Brune   PetscFunctionReturn(0);
13906a388c36SPeter Brune }
13916a388c36SPeter Brune 
13926a388c36SPeter Brune #undef __FUNCT__
1393f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchAppendOptionsPrefix"
1394e7058c64SPeter Brune /*@C
1395f1c6b773SPeter Brune    SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1396e7058c64SPeter Brune    SNES options in the database.
1397e7058c64SPeter Brune 
1398cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
1399e7058c64SPeter Brune 
1400e7058c64SPeter Brune    Input Parameters:
1401e7058c64SPeter Brune +  snes - the SNES context
1402e7058c64SPeter Brune -  prefix - the prefix to prepend to all option names
1403e7058c64SPeter Brune 
1404e7058c64SPeter Brune    Notes:
1405e7058c64SPeter Brune    A hyphen (-) must NOT be given at the beginning of the prefix name.
1406e7058c64SPeter Brune    The first character of all runtime options is AUTOMATICALLY the hyphen.
1407e7058c64SPeter Brune 
1408e7058c64SPeter Brune    Level: advanced
1409e7058c64SPeter Brune 
1410f1c6b773SPeter Brune .keywords: SNESLineSearch, append, options, prefix, database
1411e7058c64SPeter Brune 
1412e7058c64SPeter Brune .seealso: SNESGetOptionsPrefix()
1413e7058c64SPeter Brune @*/
1414f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[])
1415e7058c64SPeter Brune {
1416e7058c64SPeter Brune   PetscErrorCode ierr;
1417e7058c64SPeter Brune 
1418e7058c64SPeter Brune   PetscFunctionBegin;
1419f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1420e7058c64SPeter Brune   ierr = PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1421e7058c64SPeter Brune   PetscFunctionReturn(0);
1422e7058c64SPeter Brune }
1423e7058c64SPeter Brune 
1424e7058c64SPeter Brune #undef __FUNCT__
1425f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetOptionsPrefix"
1426e7058c64SPeter Brune /*@C
1427f1c6b773SPeter Brune    SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all
1428f1c6b773SPeter Brune    SNESLineSearch options in the database.
1429e7058c64SPeter Brune 
1430e7058c64SPeter Brune    Not Collective
1431e7058c64SPeter Brune 
1432e7058c64SPeter Brune    Input Parameter:
1433cd7522eaSPeter Brune .  linesearch - the SNESLineSearch context
1434e7058c64SPeter Brune 
1435e7058c64SPeter Brune    Output Parameter:
1436e7058c64SPeter Brune .  prefix - pointer to the prefix string used
1437e7058c64SPeter Brune 
14388e557f58SPeter Brune    Notes:
14398e557f58SPeter Brune    On the fortran side, the user should pass in a string 'prefix' of
1440e7058c64SPeter Brune    sufficient length to hold the prefix.
1441e7058c64SPeter Brune 
1442e7058c64SPeter Brune    Level: advanced
1443e7058c64SPeter Brune 
1444f1c6b773SPeter Brune .keywords: SNESLineSearch, get, options, prefix, database
1445e7058c64SPeter Brune 
1446e7058c64SPeter Brune .seealso: SNESAppendOptionsPrefix()
1447e7058c64SPeter Brune @*/
1448f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetOptionsPrefix(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 = PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1455e7058c64SPeter Brune   PetscFunctionReturn(0);
1456e7058c64SPeter Brune }
1457bf7f4e0aSPeter Brune 
1458bf7f4e0aSPeter Brune #undef __FUNCT__
1459fa0ddf94SBarry Smith #define __FUNCT__ "SNESLineSearchSetWorkVecs"
14608d359177SBarry Smith /*@C
1461fa0ddf94SBarry Smith    SNESLineSearchSetWorkVecs - Gets work vectors for the line search.
1462f40b411bSPeter Brune 
1463f40b411bSPeter Brune    Input Parameter:
1464f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
1465f40b411bSPeter Brune -  nwork - the number of work vectors
1466f40b411bSPeter Brune 
1467f40b411bSPeter Brune    Level: developer
1468f40b411bSPeter Brune 
1469fa0ddf94SBarry Smith    Developers Note: This is PETSC_EXTERN because it may be used by user written plugin SNES implementations
1470cd7522eaSPeter Brune 
1471f1c6b773SPeter Brune .keywords: SNESLineSearch, work, vector
1472f40b411bSPeter Brune 
1473fa0ddf94SBarry Smith .seealso: SNESSetWorkVecs()
1474f40b411bSPeter Brune @*/
1475fa0ddf94SBarry Smith PetscErrorCode  SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork)
1476bf7f4e0aSPeter Brune {
1477bf7f4e0aSPeter Brune   PetscErrorCode ierr;
1478bf388a1fSBarry Smith 
1479bf7f4e0aSPeter Brune   PetscFunctionBegin;
1480bf7f4e0aSPeter Brune   if (linesearch->vec_sol) {
1481bf7f4e0aSPeter Brune     ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr);
14828d359177SBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
1483bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1484bf7f4e0aSPeter Brune }
1485bf7f4e0aSPeter Brune 
1486bf7f4e0aSPeter Brune #undef __FUNCT__
1487f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetSuccess"
1488f40b411bSPeter Brune /*@
1489f1c6b773SPeter Brune    SNESLineSearchGetSuccess - Gets the success/failure status of the last line search application
1490f40b411bSPeter Brune 
1491f40b411bSPeter Brune    Input Parameters:
149278bcb3b5SPeter Brune .  linesearch - linesearch context
1493f40b411bSPeter Brune 
1494f40b411bSPeter Brune    Output Parameters:
14958e557f58SPeter Brune .  success - The success or failure status
1496f40b411bSPeter Brune 
1497cd7522eaSPeter Brune    Notes:
1498c98378a5SBarry Smith    This is typically called after SNESLineSearchApply() in order to determine if the line-search failed
1499cd7522eaSPeter Brune    (and set the SNES convergence accordingly).
1500cd7522eaSPeter Brune 
1501f40b411bSPeter Brune    Level: intermediate
1502f40b411bSPeter Brune 
1503f1c6b773SPeter Brune .seealso: SNESLineSearchSetSuccess()
1504f40b411bSPeter Brune @*/
1505f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetSuccess(SNESLineSearch linesearch, PetscBool *success)
1506bf7f4e0aSPeter Brune {
1507bf7f4e0aSPeter Brune   PetscFunctionBegin;
1508f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
15096a388c36SPeter Brune   PetscValidPointer(success, 2);
1510f5af7f23SKarl Rupp   if (success) *success = linesearch->success;
1511bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1512bf7f4e0aSPeter Brune }
1513bf7f4e0aSPeter Brune 
1514bf7f4e0aSPeter Brune #undef __FUNCT__
1515f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetSuccess"
1516f40b411bSPeter Brune /*@
1517f1c6b773SPeter Brune    SNESLineSearchSetSuccess - Sets the success/failure status of the last line search application
1518f40b411bSPeter Brune 
1519f40b411bSPeter Brune    Input Parameters:
152078bcb3b5SPeter Brune +  linesearch - linesearch context
15218e557f58SPeter Brune -  success - The success or failure status
1522f40b411bSPeter Brune 
1523cd7522eaSPeter Brune    Notes:
1524cd7522eaSPeter Brune    This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set
1525cd7522eaSPeter Brune    the success or failure of the line search method.
1526cd7522eaSPeter Brune 
152778bcb3b5SPeter Brune    Level: developer
1528f40b411bSPeter Brune 
1529f1c6b773SPeter Brune .seealso: SNESLineSearchGetSuccess()
1530f40b411bSPeter Brune @*/
1531f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetSuccess(SNESLineSearch linesearch, PetscBool success)
15326a388c36SPeter Brune {
15336a388c36SPeter Brune   PetscFunctionBegin;
15345fd66863SKarl Rupp   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
15356a388c36SPeter Brune   linesearch->success = success;
15366a388c36SPeter Brune   PetscFunctionReturn(0);
15376a388c36SPeter Brune }
15386a388c36SPeter Brune 
15396a388c36SPeter Brune #undef __FUNCT__
1540f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetVIFunctions"
15419bd66eb0SPeter Brune /*@C
1542f1c6b773SPeter Brune    SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.
15439bd66eb0SPeter Brune 
15449bd66eb0SPeter Brune    Input Parameters:
15459bd66eb0SPeter Brune +  snes - nonlinear context obtained from SNESCreate()
15469bd66eb0SPeter Brune .  projectfunc - function for projecting the function to the bounds
15479bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
15489bd66eb0SPeter Brune 
15499bd66eb0SPeter Brune    Logically Collective on SNES
15509bd66eb0SPeter Brune 
15519bd66eb0SPeter Brune    Calling sequence of projectfunc:
15529bd66eb0SPeter Brune .vb
15539bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X)
15549bd66eb0SPeter Brune .ve
15559bd66eb0SPeter Brune 
15569bd66eb0SPeter Brune     Input parameters for projectfunc:
15579bd66eb0SPeter Brune +   snes - nonlinear context
15589bd66eb0SPeter Brune -   X - current solution
15599bd66eb0SPeter Brune 
1560cd7522eaSPeter Brune     Output parameters for projectfunc:
15619bd66eb0SPeter Brune .   X - Projected solution
15629bd66eb0SPeter Brune 
15639bd66eb0SPeter Brune    Calling sequence of normfunc:
15649bd66eb0SPeter Brune .vb
15659bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm)
15669bd66eb0SPeter Brune .ve
15679bd66eb0SPeter Brune 
1568cd7522eaSPeter Brune     Input parameters for normfunc:
15699bd66eb0SPeter Brune +   snes - nonlinear context
15709bd66eb0SPeter Brune .   X - current solution
15719bd66eb0SPeter Brune -   F - current residual
15729bd66eb0SPeter Brune 
1573cd7522eaSPeter Brune     Output parameters for normfunc:
15749bd66eb0SPeter Brune .   fnorm - VI-specific norm of the function
15759bd66eb0SPeter Brune 
1576cd7522eaSPeter Brune     Notes:
1577cd7522eaSPeter Brune     The VI solvers require projection of the solution to the feasible set.  projectfunc should implement this.
1578cd7522eaSPeter Brune 
1579cd7522eaSPeter Brune     The VI solvers require special evaluation of the function norm such that the norm is only calculated
1580cd7522eaSPeter Brune     on the inactive set.  This should be implemented by normfunc.
15819bd66eb0SPeter Brune 
15829bd66eb0SPeter Brune     Level: developer
15839bd66eb0SPeter Brune 
15849bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, set, line search
15859bd66eb0SPeter Brune 
1586f1c6b773SPeter Brune .seealso: SNESLineSearchGetVIFunctions(), SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck()
15879bd66eb0SPeter Brune @*/
1588f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc)
15899bd66eb0SPeter Brune {
15909bd66eb0SPeter Brune   PetscFunctionBegin;
1591f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
15929bd66eb0SPeter Brune   if (projectfunc) linesearch->ops->viproject = projectfunc;
15939bd66eb0SPeter Brune   if (normfunc) linesearch->ops->vinorm = normfunc;
15949bd66eb0SPeter Brune   PetscFunctionReturn(0);
15959bd66eb0SPeter Brune }
15969bd66eb0SPeter Brune 
15979bd66eb0SPeter Brune #undef __FUNCT__
1598f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetVIFunctions"
15999bd66eb0SPeter Brune /*@C
1600f1c6b773SPeter Brune    SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.
16019bd66eb0SPeter Brune 
16029bd66eb0SPeter Brune    Input Parameters:
16039bd66eb0SPeter Brune .  snes - nonlinear context obtained from SNESCreate()
16049bd66eb0SPeter Brune 
16059bd66eb0SPeter Brune    Output Parameters:
16069bd66eb0SPeter Brune +  projectfunc - function for projecting the function to the bounds
16079bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
16089bd66eb0SPeter Brune 
16099bd66eb0SPeter Brune    Logically Collective on SNES
16109bd66eb0SPeter Brune 
16119bd66eb0SPeter Brune     Level: developer
16129bd66eb0SPeter Brune 
16139bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, get, line search
16149bd66eb0SPeter Brune 
1615f1c6b773SPeter Brune .seealso: SNESLineSearchSetVIFunctions(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
16169bd66eb0SPeter Brune @*/
1617f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc)
16189bd66eb0SPeter Brune {
16199bd66eb0SPeter Brune   PetscFunctionBegin;
16209bd66eb0SPeter Brune   if (projectfunc) *projectfunc = linesearch->ops->viproject;
16219bd66eb0SPeter Brune   if (normfunc) *normfunc = linesearch->ops->vinorm;
16229bd66eb0SPeter Brune   PetscFunctionReturn(0);
16239bd66eb0SPeter Brune }
16249bd66eb0SPeter Brune 
16259bd66eb0SPeter Brune #undef __FUNCT__
1626f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchRegister"
1627bf7f4e0aSPeter Brune /*@C
16281c84c290SBarry Smith   SNESLineSearchRegister - See SNESLineSearchRegister()
1629bf7f4e0aSPeter Brune 
1630bf7f4e0aSPeter Brune   Level: advanced
1631bf7f4e0aSPeter Brune @*/
1632bdf89e91SBarry Smith PetscErrorCode  SNESLineSearchRegister(const char sname[],PetscErrorCode (*function)(SNESLineSearch))
1633bf7f4e0aSPeter Brune {
1634bf7f4e0aSPeter Brune   PetscErrorCode ierr;
1635bf7f4e0aSPeter Brune 
1636bf7f4e0aSPeter Brune   PetscFunctionBegin;
1637a240a19fSJed Brown   ierr = PetscFunctionListAdd(&SNESLineSearchList,sname,function);CHKERRQ(ierr);
1638bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1639bf7f4e0aSPeter Brune }
1640