xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision ce94432eddcd14845bc7e8083b7f8ea723b9bf7d)
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 
1666b2b7091SBarry Smith .seealso:   SNESLineSearchSetPreCheck(), SNESLineSearchSetPostCheck()
1676b2b7091SBarry Smith M*/
1686b2b7091SBarry Smith 
16986d74e61SPeter Brune #undef __FUNCT__
170f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetPreCheck"
17186d74e61SPeter Brune /*@C
172f1c6b773SPeter Brune    SNESLineSearchSetPreCheck - Sets a pre-check function for the line search routine.
17386d74e61SPeter Brune 
174f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
17586d74e61SPeter Brune 
17686d74e61SPeter Brune    Input Parameters:
177f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
1786b2b7091SBarry Smith .  SNESLineSearchPreCheckFunction - [optional] function evaluation routine
17986d74e61SPeter Brune -  ctx        - [optional] user-defined context for private data for the
1800298fd71SBarry Smith                 function evaluation routine (may be NULL)
18186d74e61SPeter Brune 
18286d74e61SPeter Brune 
18386d74e61SPeter Brune    Level: intermediate
18486d74e61SPeter Brune 
18586d74e61SPeter Brune .keywords: set, linesearch, pre-check
18686d74e61SPeter Brune 
187f1c6b773SPeter Brune .seealso: SNESLineSearchSetPostCheck()
18886d74e61SPeter Brune @*/
1896b2b7091SBarry Smith PetscErrorCode  SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*SNESLineSearchPreCheckFunction)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void *ctx)
19086d74e61SPeter Brune {
1919bd66eb0SPeter Brune   PetscFunctionBegin;
192f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1936b2b7091SBarry Smith   if (SNESLineSearchPreCheckFunction) linesearch->ops->precheck = SNESLineSearchPreCheckFunction;
19486d74e61SPeter Brune   if (ctx) linesearch->precheckctx = ctx;
19586d74e61SPeter Brune   PetscFunctionReturn(0);
19686d74e61SPeter Brune }
19786d74e61SPeter Brune 
19886d74e61SPeter Brune #undef __FUNCT__
199f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetPreCheck"
20086d74e61SPeter Brune /*@C
201cd7522eaSPeter Brune    SNESLineSearchSetPreCheck - Gets the pre-check function for the line search routine.
20286d74e61SPeter Brune 
20386d74e61SPeter Brune    Input Parameters:
204f1c6b773SPeter Brune .  linesearch - the SNESLineSearch context
20586d74e61SPeter Brune 
20686d74e61SPeter Brune    Output Parameters:
20786d74e61SPeter Brune +  func       - [optional] function evaluation routine
20886d74e61SPeter Brune -  ctx        - [optional] user-defined context for private data for the
2090298fd71SBarry Smith                 function evaluation routine (may be NULL)
21086d74e61SPeter Brune 
21186d74e61SPeter Brune    Level: intermediate
21286d74e61SPeter Brune 
21386d74e61SPeter Brune .keywords: get, linesearch, pre-check
21486d74e61SPeter Brune 
215f1c6b773SPeter Brune .seealso: SNESLineSearchGetPostCheck(), SNESLineSearchSetPreCheck()
21686d74e61SPeter Brune @*/
2176b2b7091SBarry Smith PetscErrorCode  SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**SNESLineSearchPreCheckFunction)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void **ctx)
21886d74e61SPeter Brune {
2199bd66eb0SPeter Brune   PetscFunctionBegin;
220f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
2216b2b7091SBarry Smith   if (SNESLineSearchPreCheckFunction) *SNESLineSearchPreCheckFunction = linesearch->ops->precheck;
22286d74e61SPeter Brune   if (ctx) *ctx = linesearch->precheckctx;
22386d74e61SPeter Brune   PetscFunctionReturn(0);
22486d74e61SPeter Brune }
22586d74e61SPeter Brune 
2266b2b7091SBarry Smith /*MC
2276b2b7091SBarry Smith     SNESLineSearchPostheckFunction - functional form that is called after line search is complete
2286b2b7091SBarry Smith 
2296b2b7091SBarry Smith      Synopsis:
2306b2b7091SBarry Smith      #include "petscsnes.h"
2316b2b7091SBarry Smith      SNESLineSearchPostheckFunction(SNESLineSearch linesearch,Vec x,Vec y,  Vec w, *changed_y, PetscBool *changed_w);
2326b2b7091SBarry Smith 
2336b2b7091SBarry Smith      Input Parameters:
2346b2b7091SBarry Smith +      x - old solution vector
2356b2b7091SBarry Smith .      y - search direction vector
2366b2b7091SBarry Smith .      w - new solution vector
2376b2b7091SBarry Smith .      changed_y - indicates that the line search changed y
2386b2b7091SBarry Smith -      changed_w - indicates that the line search changed w
2396b2b7091SBarry Smith 
2406b2b7091SBarry Smith 
2416b2b7091SBarry Smith .seealso:   SNESLineSearchSetPreCheck(), SNESLineSearchSetPostCheck()
2426b2b7091SBarry Smith M*/
2436b2b7091SBarry Smith 
24486d74e61SPeter Brune #undef __FUNCT__
245f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetPostCheck"
24686d74e61SPeter Brune /*@C
247f1c6b773SPeter Brune    SNESLineSearchSetPostCheck - Sets a post-check function for the line search routine.
24886d74e61SPeter Brune 
249f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
25086d74e61SPeter Brune 
25186d74e61SPeter Brune    Input Parameters:
252f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
2536b2b7091SBarry Smith .  SNESLineSearchPostCheckFunction - [optional] function evaluation routine
25486d74e61SPeter Brune -  ctx        - [optional] user-defined context for private data for the
2550298fd71SBarry Smith                 function evaluation routine (may be NULL)
25686d74e61SPeter Brune 
25786d74e61SPeter Brune    Level: intermediate
25886d74e61SPeter Brune 
25986d74e61SPeter Brune .keywords: set, linesearch, post-check
26086d74e61SPeter Brune 
261f1c6b773SPeter Brune .seealso: SNESLineSearchSetPreCheck()
26286d74e61SPeter Brune @*/
2636b2b7091SBarry Smith PetscErrorCode  SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*SNESLineSearchPostCheckFunction)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void *ctx)
26486d74e61SPeter Brune {
26586d74e61SPeter Brune   PetscFunctionBegin;
266f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
2676b2b7091SBarry Smith   if (SNESLineSearchPostCheckFunction) linesearch->ops->postcheck = SNESLineSearchPostCheckFunction;
26886d74e61SPeter Brune   if (ctx) linesearch->postcheckctx = ctx;
26986d74e61SPeter Brune   PetscFunctionReturn(0);
27086d74e61SPeter Brune }
27186d74e61SPeter Brune 
27286d74e61SPeter Brune #undef __FUNCT__
273f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetPostCheck"
27486d74e61SPeter Brune /*@C
275f1c6b773SPeter Brune    SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine.
27686d74e61SPeter Brune 
27786d74e61SPeter Brune    Input Parameters:
278f1c6b773SPeter Brune .  linesearch - the SNESLineSearch context
27986d74e61SPeter Brune 
28086d74e61SPeter Brune    Output Parameters:
2816b2b7091SBarry Smith +  SNESLineSearchPostCheckFunction - [optional] function evaluation routine
28286d74e61SPeter Brune -  ctx        - [optional] user-defined context for private data for the
2830298fd71SBarry Smith                 function evaluation routine (may be NULL)
28486d74e61SPeter Brune 
28586d74e61SPeter Brune    Level: intermediate
28686d74e61SPeter Brune 
28786d74e61SPeter Brune .keywords: get, linesearch, post-check
28886d74e61SPeter Brune 
289f1c6b773SPeter Brune .seealso: SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck()
29086d74e61SPeter Brune @*/
2916b2b7091SBarry Smith PetscErrorCode  SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**SNESLineSearchPostCheckFunction)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void **ctx)
29286d74e61SPeter Brune {
2939bd66eb0SPeter Brune   PetscFunctionBegin;
294f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
2956b2b7091SBarry Smith   if (SNESLineSearchPostCheckFunction) *SNESLineSearchPostCheckFunction = linesearch->ops->postcheck;
29686d74e61SPeter Brune   if (ctx) *ctx = linesearch->postcheckctx;
29786d74e61SPeter Brune   PetscFunctionReturn(0);
29886d74e61SPeter Brune }
29986d74e61SPeter Brune 
300bf7f4e0aSPeter Brune #undef __FUNCT__
301f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPreCheck"
302f40b411bSPeter Brune /*@
303f1c6b773SPeter Brune    SNESLineSearchPreCheck - Prepares the line search for being applied.
304f40b411bSPeter Brune 
305cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
306f40b411bSPeter Brune 
307f40b411bSPeter Brune    Input Parameters:
3087b1df9c1SPeter Brune +  linesearch - The linesearch instance.
3097b1df9c1SPeter Brune .  X - The current solution
3107b1df9c1SPeter Brune -  Y - The step direction
311f40b411bSPeter Brune 
312f40b411bSPeter Brune    Output Parameters:
3138e557f58SPeter Brune .  changed - Indicator that the precheck routine has changed anything
314f40b411bSPeter Brune 
315f40b411bSPeter Brune    Level: Beginner
316f40b411bSPeter Brune 
317f1c6b773SPeter Brune .keywords: SNESLineSearch, Create
318f40b411bSPeter Brune 
319f1c6b773SPeter Brune .seealso: SNESLineSearchPostCheck()
320f40b411bSPeter Brune @*/
3217b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed)
322bf7f4e0aSPeter Brune {
323bf7f4e0aSPeter Brune   PetscErrorCode ierr;
3245fd66863SKarl Rupp 
325bf7f4e0aSPeter Brune   PetscFunctionBegin;
326bf7f4e0aSPeter Brune   *changed = PETSC_FALSE;
3276b2b7091SBarry Smith   if (linesearch->ops->precheck) {
3286b2b7091SBarry Smith     ierr = (*linesearch->ops->precheck)(linesearch, X, Y, changed, linesearch->precheckctx);CHKERRQ(ierr);
32938bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed,4);
330bf7f4e0aSPeter Brune   }
331bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
332bf7f4e0aSPeter Brune }
333bf7f4e0aSPeter Brune 
334bf7f4e0aSPeter Brune #undef __FUNCT__
335f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPostCheck"
336f40b411bSPeter Brune /*@
337f1c6b773SPeter Brune    SNESLineSearchPostCheck - Prepares the line search for being applied.
338f40b411bSPeter Brune 
339cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
340f40b411bSPeter Brune 
341f40b411bSPeter Brune    Input Parameters:
3427b1df9c1SPeter Brune +  linesearch - The linesearch context
3437b1df9c1SPeter Brune .  X - The last solution
3447b1df9c1SPeter Brune .  Y - The step direction
3457b1df9c1SPeter Brune -  W - The updated solution, W = X + lambda*Y for some lambda
346f40b411bSPeter Brune 
347f40b411bSPeter Brune    Output Parameters:
34878bcb3b5SPeter Brune +  changed_Y - Indicator if the direction Y has been changed.
34978bcb3b5SPeter Brune -  changed_W - Indicator if the new candidate solution W has been changed.
350f40b411bSPeter Brune 
351f40b411bSPeter Brune    Level: Intermediate
352f40b411bSPeter Brune 
353f1c6b773SPeter Brune .keywords: SNESLineSearch, Create
354f40b411bSPeter Brune 
355f1c6b773SPeter Brune .seealso: SNESLineSearchPreCheck()
356f40b411bSPeter Brune @*/
3577b1df9c1SPeter Brune PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch,Vec X,Vec Y,Vec W,PetscBool *changed_Y,PetscBool *changed_W)
358bf7f4e0aSPeter Brune {
359bf7f4e0aSPeter Brune   PetscErrorCode ierr;
360bf388a1fSBarry Smith 
361bf7f4e0aSPeter Brune   PetscFunctionBegin;
362bf7f4e0aSPeter Brune   *changed_Y = PETSC_FALSE;
363bf7f4e0aSPeter Brune   *changed_W = PETSC_FALSE;
3646b2b7091SBarry Smith   if (linesearch->ops->postcheck) {
3656b2b7091SBarry Smith     ierr = (*linesearch->ops->postcheck)(linesearch,X,Y,W,changed_Y,changed_W,linesearch->postcheckctx);CHKERRQ(ierr);
36638bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed_Y,5);
36738bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch,*changed_W,6);
36886d74e61SPeter Brune   }
36986d74e61SPeter Brune   PetscFunctionReturn(0);
37086d74e61SPeter Brune }
37186d74e61SPeter Brune 
37286d74e61SPeter Brune #undef __FUNCT__
373f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchPreCheckPicard"
37486d74e61SPeter Brune /*@C
37586d74e61SPeter Brune    SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration
37686d74e61SPeter Brune 
377cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
37886d74e61SPeter Brune 
37986d74e61SPeter Brune    Input Arguments:
38086d74e61SPeter Brune +  linesearch - linesearch context
38186d74e61SPeter Brune .  X - base state for this step
38286d74e61SPeter Brune .  Y - initial correction
38386d74e61SPeter Brune 
38486d74e61SPeter Brune    Output Arguments:
38586d74e61SPeter Brune +  Y - correction, possibly modified
38686d74e61SPeter Brune -  changed - flag indicating that Y was modified
38786d74e61SPeter Brune 
38886d74e61SPeter Brune    Options Database Key:
389cd7522eaSPeter Brune +  -snes_linesearch_precheck_picard - activate this routine
390cd7522eaSPeter Brune -  -snes_linesearch_precheck_picard_angle - angle
39186d74e61SPeter Brune 
39286d74e61SPeter Brune    Level: advanced
39386d74e61SPeter Brune 
39486d74e61SPeter Brune    Notes:
39586d74e61SPeter Brune    This function should be passed to SNESLineSearchSetPreCheck()
39686d74e61SPeter Brune 
39786d74e61SPeter Brune    The justification for this method involves the linear convergence of a Picard iteration
39886d74e61SPeter Brune    so the Picard linearization should be provided in place of the "Jacobian". This correction
39986d74e61SPeter Brune    is generally not useful when using a Newton linearization.
40086d74e61SPeter Brune 
40186d74e61SPeter Brune    Reference:
40286d74e61SPeter Brune    Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology.
40386d74e61SPeter Brune 
40486d74e61SPeter Brune .seealso: SNESLineSearchSetPreCheck()
40586d74e61SPeter Brune @*/
406f1c6b773SPeter Brune PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed,void *ctx)
40786d74e61SPeter Brune {
40886d74e61SPeter Brune   PetscErrorCode ierr;
40986d74e61SPeter Brune   PetscReal      angle = *(PetscReal*)linesearch->precheckctx;
41086d74e61SPeter Brune   Vec            Ylast;
41186d74e61SPeter Brune   PetscScalar    dot;
41286d74e61SPeter Brune   PetscInt       iter;
41386d74e61SPeter Brune   PetscReal      ynorm,ylastnorm,theta,angle_radians;
41486d74e61SPeter Brune   SNES           snes;
41586d74e61SPeter Brune 
41686d74e61SPeter Brune   PetscFunctionBegin;
417f1c6b773SPeter Brune   ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
41886d74e61SPeter Brune   ierr = PetscObjectQuery((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject*)&Ylast);CHKERRQ(ierr);
41986d74e61SPeter Brune   if (!Ylast) {
42086d74e61SPeter Brune     ierr = VecDuplicate(Y,&Ylast);CHKERRQ(ierr);
42186d74e61SPeter Brune     ierr = PetscObjectCompose((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject)Ylast);CHKERRQ(ierr);
42286d74e61SPeter Brune     ierr = PetscObjectDereference((PetscObject)Ylast);CHKERRQ(ierr);
42386d74e61SPeter Brune   }
42486d74e61SPeter Brune   ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr);
42586d74e61SPeter Brune   if (iter < 2) {
42686d74e61SPeter Brune     ierr     = VecCopy(Y,Ylast);CHKERRQ(ierr);
42786d74e61SPeter Brune     *changed = PETSC_FALSE;
42886d74e61SPeter Brune     PetscFunctionReturn(0);
42986d74e61SPeter Brune   }
43086d74e61SPeter Brune 
43186d74e61SPeter Brune   ierr = VecDot(Y,Ylast,&dot);CHKERRQ(ierr);
43286d74e61SPeter Brune   ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr);
43386d74e61SPeter Brune   ierr = VecNorm(Ylast,NORM_2,&ylastnorm);CHKERRQ(ierr);
43486d74e61SPeter Brune   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
43586d74e61SPeter Brune   theta         = acos((double)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0));
43686d74e61SPeter Brune   angle_radians = angle * PETSC_PI / 180.;
43786d74e61SPeter Brune   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
43886d74e61SPeter Brune     /* Modify the step Y */
43986d74e61SPeter Brune     PetscReal alpha,ydiffnorm;
44086d74e61SPeter Brune     ierr  = VecAXPY(Ylast,-1.0,Y);CHKERRQ(ierr);
44186d74e61SPeter Brune     ierr  = VecNorm(Ylast,NORM_2,&ydiffnorm);CHKERRQ(ierr);
44286d74e61SPeter Brune     alpha = ylastnorm / ydiffnorm;
44386d74e61SPeter Brune     ierr  = VecCopy(Y,Ylast);CHKERRQ(ierr);
44486d74e61SPeter Brune     ierr  = VecScale(Y,alpha);CHKERRQ(ierr);
445c69d1a72SBarry 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);
44686d74e61SPeter Brune   } else {
447c69d1a72SBarry 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);
44886d74e61SPeter Brune     ierr     = VecCopy(Y,Ylast);CHKERRQ(ierr);
44986d74e61SPeter Brune     *changed = PETSC_FALSE;
450bf7f4e0aSPeter Brune   }
451bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
452bf7f4e0aSPeter Brune }
453bf7f4e0aSPeter Brune 
454bf7f4e0aSPeter Brune #undef __FUNCT__
455f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchApply"
456f40b411bSPeter Brune /*@
457cd7522eaSPeter Brune    SNESLineSearchApply - Computes the line-search update.
458f40b411bSPeter Brune 
459f1c6b773SPeter Brune    Collective on SNESLineSearch
460f40b411bSPeter Brune 
461f40b411bSPeter Brune    Input Parameters:
4628e557f58SPeter Brune +  linesearch - The linesearch context
4638e557f58SPeter Brune .  X - The current solution
4648e557f58SPeter Brune .  F - The current function
4658e557f58SPeter Brune .  fnorm - The current norm
4668e557f58SPeter Brune -  Y - The search direction
467f40b411bSPeter Brune 
468f40b411bSPeter Brune    Output Parameters:
4698e557f58SPeter Brune +  X - The new solution
4708e557f58SPeter Brune .  F - The new function
4718e557f58SPeter Brune -  fnorm - The new function norm
472f40b411bSPeter Brune 
473cd7522eaSPeter Brune    Options Database Keys:
4743c7d6663SPeter Brune + -snes_linesearch_type - basic, bt, l2, cp, shell
475cd7522eaSPeter Brune . -snes_linesearch_monitor - Print progress of line searches
4763c7d6663SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter
477cd7522eaSPeter Brune . -snes_linesearch_norms   - Turn on/off the linesearch norms
4783c7d6663SPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess
4793c7d6663SPeter Brune - -snes_linesearch_max_it - The number of iterations for iterative line searches
480cd7522eaSPeter Brune 
481cd7522eaSPeter Brune    Notes:
482cd7522eaSPeter Brune    This is typically called from within a SNESSolve() implementation in order to
483cd7522eaSPeter Brune    help with convergence of the nonlinear method.  Various SNES types use line searches
484cd7522eaSPeter Brune    in different ways, but the overarching theme is that a line search is used to determine
485cd7522eaSPeter Brune    an optimal damping parameter of a step at each iteration of the method.  Each
486cd7522eaSPeter Brune    application of the line search may invoke SNESComputeFunction several times, and
487cd7522eaSPeter Brune    therefore may be fairly expensive.
488cd7522eaSPeter Brune 
489f40b411bSPeter Brune    Level: Intermediate
490f40b411bSPeter Brune 
491f1c6b773SPeter Brune .keywords: SNESLineSearch, Create
492f40b411bSPeter Brune 
493cd7522eaSPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchPreCheck(), SNESLineSearchPostCheck(), SNESSolve(), SNESComputeFunction()
494f40b411bSPeter Brune @*/
495bf388a1fSBarry Smith PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y)
496bf388a1fSBarry Smith {
497bf7f4e0aSPeter Brune   PetscErrorCode ierr;
498bf7f4e0aSPeter Brune 
499bf388a1fSBarry Smith   PetscFunctionBegin;
500f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
501bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
502bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
503bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
504bf7f4e0aSPeter Brune 
505bf7f4e0aSPeter Brune   linesearch->success = PETSC_TRUE;
506bf7f4e0aSPeter Brune 
507bf7f4e0aSPeter Brune   linesearch->vec_sol    = X;
508bf7f4e0aSPeter Brune   linesearch->vec_update = Y;
509bf7f4e0aSPeter Brune   linesearch->vec_func   = F;
510bf7f4e0aSPeter Brune 
511f1c6b773SPeter Brune   ierr = SNESLineSearchSetUp(linesearch);CHKERRQ(ierr);
512bf7f4e0aSPeter Brune 
513f5af7f23SKarl Rupp   if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */
514bf7f4e0aSPeter Brune 
515f5af7f23SKarl Rupp   if (fnorm) linesearch->fnorm = *fnorm;
516f5af7f23SKarl Rupp   else {
517bf7f4e0aSPeter Brune     ierr = VecNorm(F, NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
518bf7f4e0aSPeter Brune   }
519bf7f4e0aSPeter Brune 
520f1c6b773SPeter Brune   ierr = PetscLogEventBegin(SNESLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
521bf7f4e0aSPeter Brune 
522bf7f4e0aSPeter Brune   ierr = (*linesearch->ops->apply)(linesearch);CHKERRQ(ierr);
523bf7f4e0aSPeter Brune 
524f1c6b773SPeter Brune   ierr = PetscLogEventEnd(SNESLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
525bf7f4e0aSPeter Brune 
526f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
527bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
528bf7f4e0aSPeter Brune }
529bf7f4e0aSPeter Brune 
530bf7f4e0aSPeter Brune #undef __FUNCT__
531f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchDestroy"
532f40b411bSPeter Brune /*@
533f1c6b773SPeter Brune    SNESLineSearchDestroy - Destroys the line search instance.
534f40b411bSPeter Brune 
535f1c6b773SPeter Brune    Collective on SNESLineSearch
536f40b411bSPeter Brune 
537f40b411bSPeter Brune    Input Parameters:
5388e557f58SPeter Brune .  linesearch - The linesearch context
539f40b411bSPeter Brune 
540f40b411bSPeter Brune    Level: Intermediate
541f40b411bSPeter Brune 
54278bcb3b5SPeter Brune .keywords: SNESLineSearch, Destroy
543f40b411bSPeter Brune 
544cd7522eaSPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchReset(), SNESDestroy()
545f40b411bSPeter Brune @*/
546bf388a1fSBarry Smith PetscErrorCode SNESLineSearchDestroy(SNESLineSearch * linesearch)
547bf388a1fSBarry Smith {
548bf7f4e0aSPeter Brune   PetscErrorCode ierr;
549bf388a1fSBarry Smith 
550bf7f4e0aSPeter Brune   PetscFunctionBegin;
551bf7f4e0aSPeter Brune   if (!*linesearch) PetscFunctionReturn(0);
552f1c6b773SPeter Brune   PetscValidHeaderSpecific((*linesearch),SNESLINESEARCH_CLASSID,1);
553bf7f4e0aSPeter Brune   if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = 0; PetscFunctionReturn(0);}
554bf7f4e0aSPeter Brune   ierr = PetscObjectDepublish((*linesearch));CHKERRQ(ierr);
55522d28d08SBarry Smith   ierr = SNESLineSearchReset(*linesearch);CHKERRQ(ierr);
556f5af7f23SKarl Rupp   if ((*linesearch)->ops->destroy) (*linesearch)->ops->destroy(*linesearch);
557bf7f4e0aSPeter Brune   ierr = PetscViewerDestroy(&(*linesearch)->monitor);CHKERRQ(ierr);
558e7058c64SPeter Brune   ierr = PetscHeaderDestroy(linesearch);CHKERRQ(ierr);
559bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
560bf7f4e0aSPeter Brune }
561bf7f4e0aSPeter Brune 
562bf7f4e0aSPeter Brune #undef __FUNCT__
563f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetMonitor"
564f40b411bSPeter Brune /*@
565cd7522eaSPeter Brune    SNESLineSearchSetMonitor - Turns on/off printing useful information and debugging output about the line search.
566bf7f4e0aSPeter Brune 
567bf7f4e0aSPeter Brune    Input Parameters:
568bf7f4e0aSPeter Brune +  snes - nonlinear context obtained from SNESCreate()
569bf7f4e0aSPeter Brune -  flg - PETSC_TRUE to monitor the line search
570bf7f4e0aSPeter Brune 
571bf7f4e0aSPeter Brune    Logically Collective on SNES
572bf7f4e0aSPeter Brune 
573bf7f4e0aSPeter Brune    Options Database:
5748e557f58SPeter Brune .   -snes_linesearch_monitor - enables the monitor
575bf7f4e0aSPeter Brune 
576bf7f4e0aSPeter Brune    Level: intermediate
577bf7f4e0aSPeter Brune 
578bf7f4e0aSPeter Brune 
579cd7522eaSPeter Brune .seealso: SNESLineSearchGetMonitor(), PetscViewer
580bf7f4e0aSPeter Brune @*/
581f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetMonitor(SNESLineSearch linesearch, PetscBool flg)
582bf7f4e0aSPeter Brune {
583bf7f4e0aSPeter Brune   PetscErrorCode ierr;
584bf388a1fSBarry Smith 
585bf7f4e0aSPeter Brune   PetscFunctionBegin;
586bf7f4e0aSPeter Brune   if (flg && !linesearch->monitor) {
587*ce94432eSBarry Smith     ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)linesearch),"stdout",&linesearch->monitor);CHKERRQ(ierr);
588bf7f4e0aSPeter Brune   } else if (!flg && linesearch->monitor) {
589bf7f4e0aSPeter Brune     ierr = PetscViewerDestroy(&linesearch->monitor);CHKERRQ(ierr);
590bf7f4e0aSPeter Brune   }
591bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
592bf7f4e0aSPeter Brune }
593bf7f4e0aSPeter Brune 
594bf7f4e0aSPeter Brune #undef __FUNCT__
595f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetMonitor"
596f40b411bSPeter Brune /*@
597cd7522eaSPeter Brune    SNESLineSearchGetMonitor - Gets the PetscViewer instance for the line search monitor.
5986a388c36SPeter Brune 
599f40b411bSPeter Brune    Input Parameters:
6008e557f58SPeter Brune .  linesearch - linesearch context
601f40b411bSPeter Brune 
602f40b411bSPeter Brune    Input Parameters:
6038e557f58SPeter Brune .  monitor - monitor context
604f40b411bSPeter Brune 
605f40b411bSPeter Brune    Logically Collective on SNES
606f40b411bSPeter Brune 
607f40b411bSPeter Brune 
608f40b411bSPeter Brune    Options Database Keys:
6098e557f58SPeter Brune .   -snes_linesearch_monitor - enables the monitor
610f40b411bSPeter Brune 
611f40b411bSPeter Brune    Level: intermediate
612f40b411bSPeter Brune 
613f40b411bSPeter Brune 
614cd7522eaSPeter Brune .seealso: SNESLineSearchSetMonitor(), PetscViewer
615f40b411bSPeter Brune @*/
616f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetMonitor(SNESLineSearch linesearch, PetscViewer *monitor)
6176a388c36SPeter Brune {
6186a388c36SPeter Brune   PetscFunctionBegin;
619f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
6206a388c36SPeter Brune   if (monitor) {
6216a388c36SPeter Brune     PetscValidPointer(monitor, 2);
6226a388c36SPeter Brune     *monitor = linesearch->monitor;
6236a388c36SPeter Brune   }
6246a388c36SPeter Brune   PetscFunctionReturn(0);
6256a388c36SPeter Brune }
6266a388c36SPeter Brune 
6276a388c36SPeter Brune #undef __FUNCT__
628f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetFromOptions"
629f40b411bSPeter Brune /*@
630f1c6b773SPeter Brune    SNESLineSearchSetFromOptions - Sets options for the line search
631f40b411bSPeter Brune 
632f40b411bSPeter Brune    Input Parameters:
6338e557f58SPeter Brune .  linesearch - linesearch context
634f40b411bSPeter Brune 
635f40b411bSPeter Brune    Options Database Keys:
6363c7d6663SPeter Brune + -snes_linesearch_type <type> - basic, bt, l2, cp, shell
6373c7d6663SPeter Brune . -snes_linesearch_order <order> - 1, 2, 3.  Most types only support certain orders (bt supports 2 or 3)
6385a9b6599SPeter Brune . -snes_linesearch_norms   - Turn on/off the linesearch norms for the basic linesearch type
63971eef1aeSPeter Brune . -snes_linesearch_minlambda - The minimum step length
6401a9b3a06SPeter Brune . -snes_linesearch_maxstep - The maximum step size
6411a9b3a06SPeter Brune . -snes_linesearch_rtol - Relative tolerance for iterative line searches
6421a9b3a06SPeter Brune . -snes_linesearch_atol - Absolute tolerance for iterative line searches
6431a9b3a06SPeter Brune . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches
6441a9b3a06SPeter Brune . -snes_linesearch_max_it - The number of iterations for iterative line searches
645cd7522eaSPeter Brune . -snes_linesearch_monitor - Print progress of line searches
6468e557f58SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter
647cd7522eaSPeter Brune . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess.
6481a9b3a06SPeter Brune . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method
6491a9b3a06SPeter Brune - -snes_linesearch_precheck_picard_angle - Angle used in picard precheck method
650f40b411bSPeter Brune 
651f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
652f40b411bSPeter Brune 
653f40b411bSPeter Brune    Level: intermediate
654f40b411bSPeter Brune 
6553c7d6663SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetOrder(), SNESLineSearchSetType(), SNESLineSearchSetTolerances(), SNESLineSearchSetDamping(), SNESLineSearchPreCheckPicard()
656f40b411bSPeter Brune @*/
657bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch)
658bf388a1fSBarry Smith {
659bf7f4e0aSPeter Brune   PetscErrorCode ierr;
6601a4f838cSPeter Brune   const char     *deft = SNESLINESEARCHBASIC;
661bf7f4e0aSPeter Brune   char           type[256];
662bf7f4e0aSPeter Brune   PetscBool      flg, set;
663bf388a1fSBarry Smith 
664bf7f4e0aSPeter Brune   PetscFunctionBegin;
6650298fd71SBarry Smith   if (!SNESLineSearchRegisterAllCalled) {ierr = SNESLineSearchRegisterAll(NULL);CHKERRQ(ierr);}
666bf7f4e0aSPeter Brune 
667bf7f4e0aSPeter Brune   ierr = PetscObjectOptionsBegin((PetscObject)linesearch);CHKERRQ(ierr);
668f5af7f23SKarl Rupp   if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name;
6693c7d6663SPeter Brune   ierr = PetscOptionsList("-snes_linesearch_type","Linesearch type","SNESLineSearchSetType",SNESLineSearchList,deft,type,256,&flg);CHKERRQ(ierr);
670bf7f4e0aSPeter Brune   if (flg) {
671f1c6b773SPeter Brune     ierr = SNESLineSearchSetType(linesearch,type);CHKERRQ(ierr);
672bf7f4e0aSPeter Brune   } else if (!((PetscObject)linesearch)->type_name) {
673f1c6b773SPeter Brune     ierr = SNESLineSearchSetType(linesearch,deft);CHKERRQ(ierr);
674bf7f4e0aSPeter Brune   }
675bf7f4e0aSPeter Brune 
6767a35526eSPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_monitor","Print progress of line searches","SNESSNESLineSearchSetMonitor",
677bf7f4e0aSPeter Brune                           linesearch->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
678f1c6b773SPeter Brune   if (set) {ierr = SNESLineSearchSetMonitor(linesearch,flg);CHKERRQ(ierr);}
679bf7f4e0aSPeter Brune 
6801a9b3a06SPeter Brune   /* tolerances */
68171eef1aeSPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_minlambda","Minimum step length","SNESLineSearchSetTolerances",linesearch->steptol,&linesearch->steptol,0);CHKERRQ(ierr);
6821a9b3a06SPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_maxstep","Maximum step size","SNESLineSearchSetTolerances",linesearch->maxstep,&linesearch->maxstep,0);CHKERRQ(ierr);
6831a9b3a06SPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_rtol","Relative tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->rtol,&linesearch->rtol,0);CHKERRQ(ierr);
6841a9b3a06SPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_atol","Absolute tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->atol,&linesearch->atol,0);CHKERRQ(ierr);
6851a9b3a06SPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_ltol","Change in lambda tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->ltol,&linesearch->ltol,0);CHKERRQ(ierr);
6868e557f58SPeter Brune   ierr = PetscOptionsInt("-snes_linesearch_max_it","Maximum iterations for iterative line searches","SNESLineSearchSetTolerances",linesearch->max_its,&linesearch->max_its,0);CHKERRQ(ierr);
6877a35526eSPeter Brune 
6881a9b3a06SPeter Brune   /* damping parameters */
6891a9b3a06SPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_damping","Line search damping and initial step guess","SNESLineSearchSetDamping",linesearch->damping,&linesearch->damping,0);CHKERRQ(ierr);
6901a9b3a06SPeter Brune 
6911a9b3a06SPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_keeplambda","Use previous lambda as damping","SNESLineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,0);CHKERRQ(ierr);
6921a9b3a06SPeter Brune 
6931a9b3a06SPeter Brune   /* precheck */
6947a35526eSPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_precheck_picard","Use a correction that sometimes improves convergence of Picard iteration","SNESLineSearchPreCheckPicard",flg,&flg,&set);CHKERRQ(ierr);
6957a35526eSPeter Brune   if (set) {
6967a35526eSPeter Brune     if (flg) {
6977a35526eSPeter Brune       linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */
698f5af7f23SKarl Rupp 
6997a35526eSPeter Brune       ierr = PetscOptionsReal("-snes_linesearch_precheck_picard_angle","Maximum angle at which to activate the correction",
7000298fd71SBarry Smith                               "none",linesearch->precheck_picard_angle,&linesearch->precheck_picard_angle,NULL);CHKERRQ(ierr);
7017a35526eSPeter Brune       ierr = SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&linesearch->precheck_picard_angle);CHKERRQ(ierr);
7027a35526eSPeter Brune     } else {
7030298fd71SBarry Smith       ierr = SNESLineSearchSetPreCheck(linesearch,NULL,NULL);CHKERRQ(ierr);
7047a35526eSPeter Brune     }
7057a35526eSPeter Brune   }
706b000cd8dSPeter Brune   ierr = PetscOptionsInt("-snes_linesearch_order","Order of approximation used in the line search","SNESLineSearchSetOrder",linesearch->order,&linesearch->order,0);CHKERRQ(ierr);
7071a9b3a06SPeter Brune   ierr = PetscOptionsBool("-snes_linesearch_norms","Compute final norms in line search","SNESLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,0);CHKERRQ(ierr);
7087a35526eSPeter Brune 
7095a9b6599SPeter Brune   if (linesearch->ops->setfromoptions) {
7105a9b6599SPeter Brune     (*linesearch->ops->setfromoptions)(linesearch);CHKERRQ(ierr);
7115a9b6599SPeter Brune   }
7125a9b6599SPeter Brune 
713bf7f4e0aSPeter Brune   ierr = PetscObjectProcessOptionsHandlers((PetscObject)linesearch);CHKERRQ(ierr);
714bf7f4e0aSPeter Brune   ierr = PetscOptionsEnd();CHKERRQ(ierr);
715bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
716bf7f4e0aSPeter Brune }
717bf7f4e0aSPeter Brune 
718bf7f4e0aSPeter Brune #undef __FUNCT__
719f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchView"
720f40b411bSPeter Brune /*@
721cd7522eaSPeter Brune    SNESLineSearchView - Prints useful information about the line search not
722cd7522eaSPeter Brune    related to an individual call.
723f40b411bSPeter Brune 
724f40b411bSPeter Brune    Input Parameters:
7258e557f58SPeter Brune .  linesearch - linesearch context
726f40b411bSPeter Brune 
727f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
728f40b411bSPeter Brune 
729f40b411bSPeter Brune    Level: intermediate
730f40b411bSPeter Brune 
731f1c6b773SPeter Brune .seealso: SNESLineSearchCreate()
732f40b411bSPeter Brune @*/
733bf388a1fSBarry Smith PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer)
734bf388a1fSBarry Smith {
7357f1410a3SPeter Brune   PetscErrorCode ierr;
7367f1410a3SPeter Brune   PetscBool      iascii;
737bf388a1fSBarry Smith 
738bf7f4e0aSPeter Brune   PetscFunctionBegin;
7397f1410a3SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
7407f1410a3SPeter Brune   if (!viewer) {
741*ce94432eSBarry Smith     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch),&viewer);CHKERRQ(ierr);
7427f1410a3SPeter Brune   }
7437f1410a3SPeter Brune   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
7447f1410a3SPeter Brune   PetscCheckSameComm(linesearch,1,viewer,2);
745f40b411bSPeter Brune 
746251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
7477f1410a3SPeter Brune   if (iascii) {
7487f1410a3SPeter Brune     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)linesearch,viewer,"SNESLineSearch Object");CHKERRQ(ierr);
7497f1410a3SPeter Brune     if (linesearch->ops->view) {
7507f1410a3SPeter Brune       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
7517f1410a3SPeter Brune       ierr = (*linesearch->ops->view)(linesearch,viewer);CHKERRQ(ierr);
7527f1410a3SPeter Brune       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
7537f1410a3SPeter Brune     }
754c69d1a72SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep,(double)linesearch->steptol);CHKERRQ(ierr);
755c69d1a72SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol,(double)linesearch->atol,(double)linesearch->ltol);CHKERRQ(ierr);
7567f1410a3SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D\n", linesearch->max_its);CHKERRQ(ierr);
7576b2b7091SBarry Smith     if (linesearch->ops->precheck) {
7586b2b7091SBarry Smith       if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) {
7597f1410a3SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  using precheck step to speed up Picard convergence\n", linesearch->max_its);CHKERRQ(ierr);
7607f1410a3SPeter Brune       } else {
7617f1410a3SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  using user-defined precheck step\n", linesearch->max_its);CHKERRQ(ierr);
7627f1410a3SPeter Brune       }
7637f1410a3SPeter Brune     }
7646b2b7091SBarry Smith     if (linesearch->ops->postcheck) {
7657f1410a3SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer,"  using user-defined postcheck step\n", linesearch->max_its);CHKERRQ(ierr);
7667f1410a3SPeter Brune     }
7677f1410a3SPeter Brune   }
768bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
769bf7f4e0aSPeter Brune }
770bf7f4e0aSPeter Brune 
771bf7f4e0aSPeter Brune #undef __FUNCT__
772f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetType"
773ea5d4fccSPeter Brune /*@C
774f1c6b773SPeter Brune    SNESLineSearchSetType - Sets the linesearch type
775f40b411bSPeter Brune 
776f40b411bSPeter Brune    Input Parameters:
7778e557f58SPeter Brune +  linesearch - linesearch context
778f40b411bSPeter Brune -  type - The type of line search to be used
779f40b411bSPeter Brune 
780cd7522eaSPeter Brune    Available Types:
781cd7522eaSPeter Brune +  basic - Simple damping line search.
7828e557f58SPeter Brune .  bt - Backtracking line search over the L2 norm of the function
7838e557f58SPeter Brune .  l2 - Secant line search over the L2 norm of the function
7848e557f58SPeter Brune .  cp - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x)
7858e557f58SPeter Brune -  shell - User provided SNESLineSearch implementation
786cd7522eaSPeter Brune 
787f1c6b773SPeter Brune    Logically Collective on SNESLineSearch
788f40b411bSPeter Brune 
789f40b411bSPeter Brune    Level: intermediate
790f40b411bSPeter Brune 
791f40b411bSPeter Brune 
792f1c6b773SPeter Brune .seealso: SNESLineSearchCreate()
793f40b411bSPeter Brune @*/
79419fd82e9SBarry Smith PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type)
795bf7f4e0aSPeter Brune {
796f1c6b773SPeter Brune   PetscErrorCode ierr,(*r)(SNESLineSearch);
797bf7f4e0aSPeter Brune   PetscBool      match;
798bf7f4e0aSPeter Brune 
799bf7f4e0aSPeter Brune   PetscFunctionBegin;
800f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
801bf7f4e0aSPeter Brune   PetscValidCharPointer(type,2);
802bf7f4e0aSPeter Brune 
803251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)linesearch,type,&match);CHKERRQ(ierr);
804bf7f4e0aSPeter Brune   if (match) PetscFunctionReturn(0);
805bf7f4e0aSPeter Brune 
806*ce94432eSBarry Smith   ierr = PetscFunctionListFind(PetscObjectComm((PetscObject)linesearch),SNESLineSearchList,type,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr);
807bf7f4e0aSPeter Brune   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type);
808bf7f4e0aSPeter Brune   /* Destroy the previous private linesearch context */
809bf7f4e0aSPeter Brune   if (linesearch->ops->destroy) {
810bf7f4e0aSPeter Brune     ierr = (*(linesearch)->ops->destroy)(linesearch);CHKERRQ(ierr);
811f5af7f23SKarl Rupp 
8120298fd71SBarry Smith     linesearch->ops->destroy = NULL;
813bf7f4e0aSPeter Brune   }
814f1c6b773SPeter Brune   /* Reinitialize function pointers in SNESLineSearchOps structure */
815bf7f4e0aSPeter Brune   linesearch->ops->apply          = 0;
816bf7f4e0aSPeter Brune   linesearch->ops->view           = 0;
817bf7f4e0aSPeter Brune   linesearch->ops->setfromoptions = 0;
818bf7f4e0aSPeter Brune   linesearch->ops->destroy        = 0;
819bf7f4e0aSPeter Brune 
820bf7f4e0aSPeter Brune   ierr = PetscObjectChangeTypeName((PetscObject)linesearch,type);CHKERRQ(ierr);
821bf7f4e0aSPeter Brune   ierr = (*r)(linesearch);CHKERRQ(ierr);
822bf7f4e0aSPeter Brune #if defined(PETSC_HAVE_AMS)
823bf7f4e0aSPeter Brune   if (PetscAMSPublishAll) {
824bf7f4e0aSPeter Brune     ierr = PetscObjectAMSPublish((PetscObject)linesearch);CHKERRQ(ierr);
825bf7f4e0aSPeter Brune   }
826bf7f4e0aSPeter Brune #endif
827bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
828bf7f4e0aSPeter Brune }
829bf7f4e0aSPeter Brune 
830bf7f4e0aSPeter Brune #undef __FUNCT__
831f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetSNES"
832f40b411bSPeter Brune /*@
83378bcb3b5SPeter Brune    SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation.
834f40b411bSPeter Brune 
835f40b411bSPeter Brune    Input Parameters:
8368e557f58SPeter Brune +  linesearch - linesearch context
837f40b411bSPeter Brune -  snes - The snes instance
838f40b411bSPeter Brune 
83978bcb3b5SPeter Brune    Level: developer
84078bcb3b5SPeter Brune 
84178bcb3b5SPeter Brune    Notes:
84278bcb3b5SPeter Brune    This happens automatically when the line search is gotten/created with
84378bcb3b5SPeter Brune    SNESGetSNESLineSearch().  This routine is therefore mainly called within SNES
84478bcb3b5SPeter Brune    implementations.
845f40b411bSPeter Brune 
8468141a3b9SPeter Brune    Level: developer
847f40b411bSPeter Brune 
848cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
849f40b411bSPeter Brune @*/
850bf388a1fSBarry Smith PetscErrorCode  SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes)
851bf388a1fSBarry Smith {
852bf7f4e0aSPeter Brune   PetscFunctionBegin;
853f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
854bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,2);
855bf7f4e0aSPeter Brune   linesearch->snes = snes;
856bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
857bf7f4e0aSPeter Brune }
858bf7f4e0aSPeter Brune 
859bf7f4e0aSPeter Brune #undef __FUNCT__
860f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetSNES"
861f40b411bSPeter Brune /*@
8628141a3b9SPeter Brune    SNESLineSearchGetSNES - Gets the SNES instance associated with the line search.
8638141a3b9SPeter Brune    Having an associated SNES is necessary because most line search implementations must be able to
8648141a3b9SPeter Brune    evaluate the function using SNESComputeFunction() for the associated SNES.  This routine
8658141a3b9SPeter Brune    is used in the line search implementations when one must get this associated SNES instance.
866f40b411bSPeter Brune 
867f40b411bSPeter Brune    Input Parameters:
8688e557f58SPeter Brune .  linesearch - linesearch context
869f40b411bSPeter Brune 
870f40b411bSPeter Brune    Output Parameters:
871f40b411bSPeter Brune .  snes - The snes instance
872f40b411bSPeter Brune 
8738141a3b9SPeter Brune    Level: developer
874f40b411bSPeter Brune 
875cd7522eaSPeter Brune .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
876f40b411bSPeter Brune @*/
877bf388a1fSBarry Smith PetscErrorCode  SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes)
878bf388a1fSBarry Smith {
879bf7f4e0aSPeter Brune   PetscFunctionBegin;
880f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
8816a388c36SPeter Brune   PetscValidPointer(snes, 2);
882bf7f4e0aSPeter Brune   *snes = linesearch->snes;
883bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
884bf7f4e0aSPeter Brune }
885bf7f4e0aSPeter Brune 
8866a388c36SPeter Brune #undef __FUNCT__
887f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetLambda"
888f40b411bSPeter Brune /*@
889f1c6b773SPeter Brune    SNESLineSearchGetLambda - Gets the last linesearch steplength discovered.
890f40b411bSPeter Brune 
891f40b411bSPeter Brune    Input Parameters:
8928e557f58SPeter Brune .  linesearch - linesearch context
893f40b411bSPeter Brune 
894f40b411bSPeter Brune    Output Parameters:
895cd7522eaSPeter Brune .  lambda - The last steplength computed during SNESLineSearchApply()
896f40b411bSPeter Brune 
89778bcb3b5SPeter Brune    Level: advanced
89878bcb3b5SPeter Brune 
8998e557f58SPeter Brune    Notes:
9008e557f58SPeter Brune    This is useful in methods where the solver is ill-scaled and
90178bcb3b5SPeter Brune    requires some adaptive notion of the difference in scale between the
90278bcb3b5SPeter Brune    solution and the function.  For instance, SNESQN may be scaled by the
90378bcb3b5SPeter Brune    line search lambda using the argument -snes_qn_scaling ls.
90478bcb3b5SPeter Brune 
905f40b411bSPeter Brune 
906cd7522eaSPeter Brune .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetDamping(), SNESLineSearchApply()
907f40b411bSPeter Brune @*/
908f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda)
9096a388c36SPeter Brune {
9106a388c36SPeter Brune   PetscFunctionBegin;
911f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
9126a388c36SPeter Brune   PetscValidPointer(lambda, 2);
9136a388c36SPeter Brune   *lambda = linesearch->lambda;
9146a388c36SPeter Brune   PetscFunctionReturn(0);
9156a388c36SPeter Brune }
9166a388c36SPeter Brune 
9176a388c36SPeter Brune #undef __FUNCT__
918f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetLambda"
919f40b411bSPeter Brune /*@
920f1c6b773SPeter Brune    SNESLineSearchSetLambda - Sets the linesearch steplength.
921f40b411bSPeter Brune 
922f40b411bSPeter Brune    Input Parameters:
9238e557f58SPeter Brune +  linesearch - linesearch context
924f40b411bSPeter Brune -  lambda - The last steplength.
925f40b411bSPeter Brune 
926cd7522eaSPeter Brune    Notes:
927cd7522eaSPeter Brune    This routine is typically used within implementations of SNESLineSearchApply
928cd7522eaSPeter Brune    to set the final steplength.  This routine (and SNESLineSearchGetLambda()) were
929cd7522eaSPeter Brune    added in order to facilitate Quasi-Newton methods that use the previous steplength
930cd7522eaSPeter Brune    as an inner scaling parameter.
931cd7522eaSPeter Brune 
93278bcb3b5SPeter Brune    Level: advanced
933f40b411bSPeter Brune 
934f1c6b773SPeter Brune .seealso: SNESLineSearchGetLambda()
935f40b411bSPeter Brune @*/
936f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
9376a388c36SPeter Brune {
9386a388c36SPeter Brune   PetscFunctionBegin;
939f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
9406a388c36SPeter Brune   linesearch->lambda = lambda;
9416a388c36SPeter Brune   PetscFunctionReturn(0);
9426a388c36SPeter Brune }
9436a388c36SPeter Brune 
9446a388c36SPeter Brune #undef  __FUNCT__
945f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetTolerances"
946f40b411bSPeter Brune /*@
9473c7d6663SPeter Brune    SNESLineSearchGetTolerances - Gets the tolerances for the linesearch.  These include
94878bcb3b5SPeter Brune    tolerances for the relative and absolute change in the function norm, the change
94978bcb3b5SPeter Brune    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
95078bcb3b5SPeter Brune    and the maximum number of iterations the line search procedure may take.
951f40b411bSPeter Brune 
952f40b411bSPeter Brune    Input Parameters:
9538e557f58SPeter Brune .  linesearch - linesearch context
954f40b411bSPeter Brune 
955f40b411bSPeter Brune    Output Parameters:
956516fe3c3SPeter Brune +  steptol - The minimum steplength
9576cc8e53bSPeter Brune .  maxstep - The maximum steplength
958516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
959516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
960516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
961516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
962f40b411bSPeter Brune 
96378bcb3b5SPeter Brune    Level: intermediate
96478bcb3b5SPeter Brune 
96578bcb3b5SPeter Brune    Notes:
96678bcb3b5SPeter Brune    Different line searches may implement these parameters slightly differently as
9673c7d6663SPeter Brune    the type requires.
968516fe3c3SPeter Brune 
969f1c6b773SPeter Brune .seealso: SNESLineSearchSetTolerances()
970f40b411bSPeter Brune @*/
971f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
9726a388c36SPeter Brune {
9736a388c36SPeter Brune   PetscFunctionBegin;
974f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
975516fe3c3SPeter Brune   if (steptol) {
9766a388c36SPeter Brune     PetscValidPointer(steptol, 2);
9776a388c36SPeter Brune     *steptol = linesearch->steptol;
978516fe3c3SPeter Brune   }
979516fe3c3SPeter Brune   if (maxstep) {
980516fe3c3SPeter Brune     PetscValidPointer(maxstep, 3);
981516fe3c3SPeter Brune     *maxstep = linesearch->maxstep;
982516fe3c3SPeter Brune   }
983516fe3c3SPeter Brune   if (rtol) {
984516fe3c3SPeter Brune     PetscValidPointer(rtol, 4);
985516fe3c3SPeter Brune     *rtol = linesearch->rtol;
986516fe3c3SPeter Brune   }
987516fe3c3SPeter Brune   if (atol) {
988516fe3c3SPeter Brune     PetscValidPointer(atol, 5);
989516fe3c3SPeter Brune     *atol = linesearch->atol;
990516fe3c3SPeter Brune   }
991516fe3c3SPeter Brune   if (ltol) {
992516fe3c3SPeter Brune     PetscValidPointer(ltol, 6);
993516fe3c3SPeter Brune     *ltol = linesearch->ltol;
994516fe3c3SPeter Brune   }
995516fe3c3SPeter Brune   if (max_its) {
996516fe3c3SPeter Brune     PetscValidPointer(max_its, 7);
997516fe3c3SPeter Brune     *max_its = linesearch->max_its;
998516fe3c3SPeter Brune   }
9996a388c36SPeter Brune   PetscFunctionReturn(0);
10006a388c36SPeter Brune }
10016a388c36SPeter Brune 
10026a388c36SPeter Brune #undef  __FUNCT__
1003f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetTolerances"
1004f40b411bSPeter Brune /*@
10053c7d6663SPeter Brune    SNESLineSearchSetTolerances -  Gets the tolerances for the linesearch.  These include
100678bcb3b5SPeter Brune    tolerances for the relative and absolute change in the function norm, the change
100778bcb3b5SPeter Brune    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
100878bcb3b5SPeter Brune    and the maximum number of iterations the line search procedure may take.
1009f40b411bSPeter Brune 
1010f40b411bSPeter Brune    Input Parameters:
10118e557f58SPeter Brune +  linesearch - linesearch context
1012516fe3c3SPeter Brune .  steptol - The minimum steplength
10136cc8e53bSPeter Brune .  maxstep - The maximum steplength
1014516fe3c3SPeter Brune .  rtol    - The relative tolerance for iterative line searches
1015516fe3c3SPeter Brune .  atol    - The absolute tolerance for iterative line searches
1016516fe3c3SPeter Brune .  ltol    - The change in lambda tolerance for iterative line searches
1017516fe3c3SPeter Brune -  max_it  - The maximum number of iterations of the line search
1018f40b411bSPeter Brune 
101978bcb3b5SPeter Brune    Notes:
10203c7d6663SPeter Brune    The user may choose to not set any of the tolerances using PETSC_DEFAULT in
102178bcb3b5SPeter Brune    place of an argument.
1022f40b411bSPeter Brune 
102378bcb3b5SPeter Brune    Level: intermediate
1024516fe3c3SPeter Brune 
1025f1c6b773SPeter Brune .seealso: SNESLineSearchGetTolerances()
1026f40b411bSPeter Brune @*/
1027f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its)
10286a388c36SPeter Brune {
10296a388c36SPeter Brune   PetscFunctionBegin;
1030f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1031d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,steptol,2);
1032d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,maxstep,3);
1033d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,rtol,4);
1034d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,atol,5);
1035d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch,ltol,6);
1036d3952184SSatish Balay   PetscValidLogicalCollectiveInt(linesearch,max_its,7);
1037d3952184SSatish Balay 
1038d3952184SSatish Balay   if (steptol!= PETSC_DEFAULT) {
1039*ce94432eSBarry Smith     if (steptol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Minimum step length %14.12e must be non-negative",(double)steptol);
10406a388c36SPeter Brune     linesearch->steptol = steptol;
1041d3952184SSatish Balay   }
1042d3952184SSatish Balay 
1043d3952184SSatish Balay   if (maxstep!= PETSC_DEFAULT) {
1044*ce94432eSBarry Smith     if (maxstep < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum step length %14.12e must be non-negative",(double)maxstep);
1045516fe3c3SPeter Brune     linesearch->maxstep = maxstep;
1046d3952184SSatish Balay   }
1047d3952184SSatish Balay 
1048d3952184SSatish Balay   if (rtol != PETSC_DEFAULT) {
1049*ce94432eSBarry 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);
1050516fe3c3SPeter Brune     linesearch->rtol = rtol;
1051d3952184SSatish Balay   }
1052d3952184SSatish Balay 
1053d3952184SSatish Balay   if (atol != PETSC_DEFAULT) {
1054*ce94432eSBarry Smith     if (atol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %14.12e must be non-negative",(double)atol);
1055516fe3c3SPeter Brune     linesearch->atol = atol;
1056d3952184SSatish Balay   }
1057d3952184SSatish Balay 
1058d3952184SSatish Balay   if (ltol != PETSC_DEFAULT) {
1059*ce94432eSBarry Smith     if (ltol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Labmda tolerance %14.12e must be non-negative",(double)ltol);
1060516fe3c3SPeter Brune     linesearch->ltol = ltol;
1061d3952184SSatish Balay   }
1062d3952184SSatish Balay 
1063d3952184SSatish Balay   if (max_its != PETSC_DEFAULT) {
1064*ce94432eSBarry Smith     if (max_its < 0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",max_its);
1065516fe3c3SPeter Brune     linesearch->max_its = max_its;
1066d3952184SSatish Balay   }
10676a388c36SPeter Brune   PetscFunctionReturn(0);
10686a388c36SPeter Brune }
10696a388c36SPeter Brune 
10706a388c36SPeter Brune #undef __FUNCT__
1071f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetDamping"
1072f40b411bSPeter Brune /*@
1073f1c6b773SPeter Brune    SNESLineSearchGetDamping - Gets the line search damping parameter.
1074f40b411bSPeter Brune 
1075f40b411bSPeter Brune    Input Parameters:
10768e557f58SPeter Brune .  linesearch - linesearch context
1077f40b411bSPeter Brune 
1078f40b411bSPeter Brune    Output Parameters:
10798e557f58SPeter Brune .  damping - The damping parameter
1080f40b411bSPeter Brune 
108178bcb3b5SPeter Brune    Level: advanced
1082f40b411bSPeter Brune 
108378bcb3b5SPeter Brune .seealso: SNESLineSearchGetStepTolerance(), SNESQN
1084f40b411bSPeter Brune @*/
1085f40b411bSPeter Brune 
1086f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping)
10876a388c36SPeter Brune {
10886a388c36SPeter Brune   PetscFunctionBegin;
1089f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
10906a388c36SPeter Brune   PetscValidPointer(damping, 2);
10916a388c36SPeter Brune   *damping = linesearch->damping;
10926a388c36SPeter Brune   PetscFunctionReturn(0);
10936a388c36SPeter Brune }
10946a388c36SPeter Brune 
10956a388c36SPeter Brune #undef __FUNCT__
1096f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetDamping"
1097f40b411bSPeter Brune /*@
1098f1c6b773SPeter Brune    SNESLineSearchSetDamping - Sets the line search damping paramter.
1099f40b411bSPeter Brune 
1100f40b411bSPeter Brune    Input Parameters:
110178bcb3b5SPeter Brune .  linesearch - linesearch context
110278bcb3b5SPeter Brune .  damping - The damping parameter
1103f40b411bSPeter Brune 
1104f40b411bSPeter Brune    Level: intermediate
1105f40b411bSPeter Brune 
1106cd7522eaSPeter Brune    Notes:
1107cd7522eaSPeter Brune    The basic line search merely takes the update step scaled by the damping parameter.
1108cd7522eaSPeter Brune    The use of the damping parameter in the l2 and cp line searches is much more subtle;
110978bcb3b5SPeter Brune    it is used as a starting point in calculating the secant step. However, the eventual
1110cd7522eaSPeter Brune    step may be of greater length than the damping parameter.  In the bt line search it is
1111cd7522eaSPeter Brune    used as the maximum possible step length, as the bt line search only backtracks.
1112cd7522eaSPeter Brune 
1113f1c6b773SPeter Brune .seealso: SNESLineSearchGetDamping()
1114f40b411bSPeter Brune @*/
1115f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping)
11166a388c36SPeter Brune {
11176a388c36SPeter Brune   PetscFunctionBegin;
1118f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
11196a388c36SPeter Brune   linesearch->damping = damping;
11206a388c36SPeter Brune   PetscFunctionReturn(0);
11216a388c36SPeter Brune }
11226a388c36SPeter Brune 
11236a388c36SPeter Brune #undef __FUNCT__
112459405d5eSPeter Brune #define __FUNCT__ "SNESLineSearchGetOrder"
112559405d5eSPeter Brune /*@
112659405d5eSPeter Brune    SNESLineSearchGetOrder - Gets the line search approximation order.
112759405d5eSPeter Brune 
112859405d5eSPeter Brune    Input Parameters:
112978bcb3b5SPeter Brune .  linesearch - linesearch context
113059405d5eSPeter Brune 
113159405d5eSPeter Brune    Output Parameters:
11328e557f58SPeter Brune .  order - The order
113359405d5eSPeter Brune 
113478bcb3b5SPeter Brune    Possible Values for order:
11353c7d6663SPeter Brune +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
11363c7d6663SPeter Brune .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
11373c7d6663SPeter Brune -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
113878bcb3b5SPeter Brune 
113959405d5eSPeter Brune    Level: intermediate
114059405d5eSPeter Brune 
114159405d5eSPeter Brune .seealso: SNESLineSearchSetOrder()
114259405d5eSPeter Brune @*/
114359405d5eSPeter Brune 
1144b000cd8dSPeter Brune PetscErrorCode  SNESLineSearchGetOrder(SNESLineSearch linesearch,PetscInt *order)
114559405d5eSPeter Brune {
114659405d5eSPeter Brune   PetscFunctionBegin;
114759405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
114859405d5eSPeter Brune   PetscValidPointer(order, 2);
114959405d5eSPeter Brune   *order = linesearch->order;
115059405d5eSPeter Brune   PetscFunctionReturn(0);
115159405d5eSPeter Brune }
115259405d5eSPeter Brune 
115359405d5eSPeter Brune #undef __FUNCT__
115459405d5eSPeter Brune #define __FUNCT__ "SNESLineSearchSetOrder"
115559405d5eSPeter Brune /*@
115659405d5eSPeter Brune    SNESLineSearchSetOrder - Sets the line search damping paramter.
115759405d5eSPeter Brune 
115859405d5eSPeter Brune    Input Parameters:
115978bcb3b5SPeter Brune .  linesearch - linesearch context
116078bcb3b5SPeter Brune .  order - The damping parameter
116159405d5eSPeter Brune 
116259405d5eSPeter Brune    Level: intermediate
116359405d5eSPeter Brune 
116478bcb3b5SPeter Brune    Possible Values for order:
11653c7d6663SPeter Brune +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
11663c7d6663SPeter Brune .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
11673c7d6663SPeter Brune -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
116878bcb3b5SPeter Brune 
116959405d5eSPeter Brune    Notes:
117059405d5eSPeter Brune    Variable orders are supported by the following line searches:
117178bcb3b5SPeter Brune +  bt - cubic and quadratic
117278bcb3b5SPeter Brune -  cp - linear and quadratic
117359405d5eSPeter Brune 
117459405d5eSPeter Brune .seealso: SNESLineSearchGetOrder()
117559405d5eSPeter Brune @*/
1176b000cd8dSPeter Brune PetscErrorCode  SNESLineSearchSetOrder(SNESLineSearch linesearch,PetscInt order)
117759405d5eSPeter Brune {
117859405d5eSPeter Brune   PetscFunctionBegin;
117959405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
118059405d5eSPeter Brune   linesearch->order = order;
118159405d5eSPeter Brune   PetscFunctionReturn(0);
118259405d5eSPeter Brune }
118359405d5eSPeter Brune 
118459405d5eSPeter Brune #undef __FUNCT__
1185f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetNorms"
1186f40b411bSPeter Brune /*@
1187f1c6b773SPeter Brune    SNESLineSearchGetNorms - Gets the norms for for X, Y, and F.
1188f40b411bSPeter Brune 
1189f40b411bSPeter Brune    Input Parameters:
119078bcb3b5SPeter Brune .  linesearch - linesearch context
1191f40b411bSPeter Brune 
1192f40b411bSPeter Brune    Output Parameters:
1193f40b411bSPeter Brune +  xnorm - The norm of the current solution
1194f40b411bSPeter Brune .  fnorm - The norm of the current function
1195f40b411bSPeter Brune -  ynorm - The norm of the current update
1196f40b411bSPeter Brune 
1197cd7522eaSPeter Brune    Notes:
1198cd7522eaSPeter Brune    This function is mainly called from SNES implementations.
1199cd7522eaSPeter Brune 
120078bcb3b5SPeter Brune    Level: developer
1201f40b411bSPeter Brune 
1202f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms() SNESLineSearchGetVecs()
1203f40b411bSPeter Brune @*/
1204f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm)
1205bf7f4e0aSPeter Brune {
1206bf7f4e0aSPeter Brune   PetscFunctionBegin;
1207f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1208f5af7f23SKarl Rupp   if (xnorm) *xnorm = linesearch->xnorm;
1209f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
1210f5af7f23SKarl Rupp   if (ynorm) *ynorm = linesearch->ynorm;
1211bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1212bf7f4e0aSPeter Brune }
1213bf7f4e0aSPeter Brune 
1214e7058c64SPeter Brune #undef __FUNCT__
1215f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetNorms"
1216f40b411bSPeter Brune /*@
1217f1c6b773SPeter Brune    SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F.
1218f40b411bSPeter Brune 
1219f40b411bSPeter Brune    Input Parameters:
122078bcb3b5SPeter Brune +  linesearch - linesearch context
1221f40b411bSPeter Brune .  xnorm - The norm of the current solution
1222f40b411bSPeter Brune .  fnorm - The norm of the current function
1223f40b411bSPeter Brune -  ynorm - The norm of the current update
1224f40b411bSPeter Brune 
122578bcb3b5SPeter Brune    Level: advanced
1226f40b411bSPeter Brune 
1227f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1228f40b411bSPeter Brune @*/
1229f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
12306a388c36SPeter Brune {
12316a388c36SPeter Brune   PetscFunctionBegin;
1232f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
12336a388c36SPeter Brune   linesearch->xnorm = xnorm;
12346a388c36SPeter Brune   linesearch->fnorm = fnorm;
12356a388c36SPeter Brune   linesearch->ynorm = ynorm;
12366a388c36SPeter Brune   PetscFunctionReturn(0);
12376a388c36SPeter Brune }
12386a388c36SPeter Brune 
12396a388c36SPeter Brune #undef __FUNCT__
1240f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchComputeNorms"
1241f40b411bSPeter Brune /*@
1242f1c6b773SPeter Brune    SNESLineSearchComputeNorms - Computes the norms of X, F, and Y.
1243f40b411bSPeter Brune 
1244f40b411bSPeter Brune    Input Parameters:
124578bcb3b5SPeter Brune .  linesearch - linesearch context
1246f40b411bSPeter Brune 
1247f40b411bSPeter Brune    Options Database Keys:
12488e557f58SPeter Brune .   -snes_linesearch_norms - turn norm computation on or off
1249f40b411bSPeter Brune 
1250f40b411bSPeter Brune    Level: intermediate
1251f40b411bSPeter Brune 
125278bcb3b5SPeter Brune .seealso: SNESLineSearchGetNorms, SNESLineSearchSetNorms(), SNESLineSearchSetComputeNorms()
1253f40b411bSPeter Brune @*/
1254f1c6b773SPeter Brune PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
12556a388c36SPeter Brune {
12566a388c36SPeter Brune   PetscErrorCode ierr;
12579bd66eb0SPeter Brune   SNES           snes;
1258bf388a1fSBarry Smith 
12596a388c36SPeter Brune   PetscFunctionBegin;
12606a388c36SPeter Brune   if (linesearch->norms) {
12619bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
1262f1c6b773SPeter Brune       ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
12639bd66eb0SPeter Brune       ierr = VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
12649bd66eb0SPeter Brune       ierr = VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
12659bd66eb0SPeter Brune       ierr = (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);CHKERRQ(ierr);
12669bd66eb0SPeter Brune     } else {
12676a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_func,   NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
12686a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_sol,    NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
12696a388c36SPeter Brune       ierr = VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
12706a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_func,     NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
12716a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_sol,      NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
12726a388c36SPeter Brune       ierr = VecNormEnd(linesearch->vec_update,   NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
12736a388c36SPeter Brune     }
12749bd66eb0SPeter Brune   }
12756a388c36SPeter Brune   PetscFunctionReturn(0);
12766a388c36SPeter Brune }
12776a388c36SPeter Brune 
12786f263ca3SPeter Brune #undef __FUNCT__
12796f263ca3SPeter Brune #define __FUNCT__ "SNESLineSearchSetComputeNorms"
12806f263ca3SPeter Brune /*@
12816f263ca3SPeter Brune    SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search.
12826f263ca3SPeter Brune 
12836f263ca3SPeter Brune    Input Parameters:
128478bcb3b5SPeter Brune +  linesearch  - linesearch context
128578bcb3b5SPeter Brune -  flg  - indicates whether or not to compute norms
12866f263ca3SPeter Brune 
12876f263ca3SPeter Brune    Options Database Keys:
12888e557f58SPeter Brune .   -snes_linesearch_norms - turn norm computation on or off
12896f263ca3SPeter Brune 
12906f263ca3SPeter Brune    Notes:
12911a4f838cSPeter Brune    This is most relevant to the SNESLINESEARCHBASIC line search type.
12926f263ca3SPeter Brune 
12936f263ca3SPeter Brune    Level: intermediate
12946f263ca3SPeter Brune 
12951a4f838cSPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetNorms(), SNESLineSearchComputeNorms(), SNESLINESEARCHBASIC
12966f263ca3SPeter Brune @*/
12976f263ca3SPeter Brune PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg)
12986f263ca3SPeter Brune {
12996f263ca3SPeter Brune   PetscFunctionBegin;
13006f263ca3SPeter Brune   linesearch->norms = flg;
13016f263ca3SPeter Brune   PetscFunctionReturn(0);
13026f263ca3SPeter Brune }
13036f263ca3SPeter Brune 
13046a388c36SPeter Brune #undef __FUNCT__
1305f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetVecs"
1306f40b411bSPeter Brune /*@
1307f1c6b773SPeter Brune    SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context
1308f40b411bSPeter Brune 
1309f40b411bSPeter Brune    Input Parameters:
131078bcb3b5SPeter Brune .  linesearch - linesearch context
1311f40b411bSPeter Brune 
1312f40b411bSPeter Brune    Output Parameters:
1313f40b411bSPeter Brune +  X - The old solution
1314f40b411bSPeter Brune .  F - The old function
1315f40b411bSPeter Brune .  Y - The search direction
1316f40b411bSPeter Brune .  W - The new solution
1317f40b411bSPeter Brune -  G - The new function
1318f40b411bSPeter Brune 
131978bcb3b5SPeter Brune    Level: advanced
1320f40b411bSPeter Brune 
1321f1c6b773SPeter Brune .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1322f40b411bSPeter Brune @*/
1323bf388a1fSBarry Smith PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G)
1324bf388a1fSBarry Smith {
13256a388c36SPeter Brune   PetscFunctionBegin;
1326f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
13276a388c36SPeter Brune   if (X) {
13286a388c36SPeter Brune     PetscValidPointer(X, 2);
13296a388c36SPeter Brune     *X = linesearch->vec_sol;
13306a388c36SPeter Brune   }
13316a388c36SPeter Brune   if (F) {
13326a388c36SPeter Brune     PetscValidPointer(F, 3);
13336a388c36SPeter Brune     *F = linesearch->vec_func;
13346a388c36SPeter Brune   }
13356a388c36SPeter Brune   if (Y) {
13366a388c36SPeter Brune     PetscValidPointer(Y, 4);
13376a388c36SPeter Brune     *Y = linesearch->vec_update;
13386a388c36SPeter Brune   }
13396a388c36SPeter Brune   if (W) {
13406a388c36SPeter Brune     PetscValidPointer(W, 5);
13416a388c36SPeter Brune     *W = linesearch->vec_sol_new;
13426a388c36SPeter Brune   }
13436a388c36SPeter Brune   if (G) {
13446a388c36SPeter Brune     PetscValidPointer(G, 6);
13456a388c36SPeter Brune     *G = linesearch->vec_func_new;
13466a388c36SPeter Brune   }
13476a388c36SPeter Brune   PetscFunctionReturn(0);
13486a388c36SPeter Brune }
13496a388c36SPeter Brune 
13506a388c36SPeter Brune #undef __FUNCT__
1351f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetVecs"
1352f40b411bSPeter Brune /*@
1353f1c6b773SPeter Brune    SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context
1354f40b411bSPeter Brune 
1355f40b411bSPeter Brune    Input Parameters:
135678bcb3b5SPeter Brune +  linesearch - linesearch context
1357f40b411bSPeter Brune .  X - The old solution
1358f40b411bSPeter Brune .  F - The old function
1359f40b411bSPeter Brune .  Y - The search direction
1360f40b411bSPeter Brune .  W - The new solution
1361f40b411bSPeter Brune -  G - The new function
1362f40b411bSPeter Brune 
136378bcb3b5SPeter Brune    Level: advanced
1364f40b411bSPeter Brune 
1365f1c6b773SPeter Brune .seealso: SNESLineSearchSetNorms(), SNESLineSearchGetVecs()
1366f40b411bSPeter Brune @*/
1367bf388a1fSBarry Smith PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G)
1368bf388a1fSBarry Smith {
13696a388c36SPeter Brune   PetscFunctionBegin;
1370f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
13716a388c36SPeter Brune   if (X) {
13726a388c36SPeter Brune     PetscValidHeaderSpecific(X,VEC_CLASSID,2);
13736a388c36SPeter Brune     linesearch->vec_sol = X;
13746a388c36SPeter Brune   }
13756a388c36SPeter Brune   if (F) {
13766a388c36SPeter Brune     PetscValidHeaderSpecific(F,VEC_CLASSID,3);
13776a388c36SPeter Brune     linesearch->vec_func = F;
13786a388c36SPeter Brune   }
13796a388c36SPeter Brune   if (Y) {
13806a388c36SPeter Brune     PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
13816a388c36SPeter Brune     linesearch->vec_update = Y;
13826a388c36SPeter Brune   }
13836a388c36SPeter Brune   if (W) {
13846a388c36SPeter Brune     PetscValidHeaderSpecific(W,VEC_CLASSID,5);
13856a388c36SPeter Brune     linesearch->vec_sol_new = W;
13866a388c36SPeter Brune   }
13876a388c36SPeter Brune   if (G) {
13886a388c36SPeter Brune     PetscValidHeaderSpecific(G,VEC_CLASSID,6);
13896a388c36SPeter Brune     linesearch->vec_func_new = G;
13906a388c36SPeter Brune   }
13916a388c36SPeter Brune   PetscFunctionReturn(0);
13926a388c36SPeter Brune }
13936a388c36SPeter Brune 
13946a388c36SPeter Brune #undef __FUNCT__
1395f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchAppendOptionsPrefix"
1396e7058c64SPeter Brune /*@C
1397f1c6b773SPeter Brune    SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1398e7058c64SPeter Brune    SNES options in the database.
1399e7058c64SPeter Brune 
1400cd7522eaSPeter Brune    Logically Collective on SNESLineSearch
1401e7058c64SPeter Brune 
1402e7058c64SPeter Brune    Input Parameters:
1403e7058c64SPeter Brune +  snes - the SNES context
1404e7058c64SPeter Brune -  prefix - the prefix to prepend to all option names
1405e7058c64SPeter Brune 
1406e7058c64SPeter Brune    Notes:
1407e7058c64SPeter Brune    A hyphen (-) must NOT be given at the beginning of the prefix name.
1408e7058c64SPeter Brune    The first character of all runtime options is AUTOMATICALLY the hyphen.
1409e7058c64SPeter Brune 
1410e7058c64SPeter Brune    Level: advanced
1411e7058c64SPeter Brune 
1412f1c6b773SPeter Brune .keywords: SNESLineSearch, append, options, prefix, database
1413e7058c64SPeter Brune 
1414e7058c64SPeter Brune .seealso: SNESGetOptionsPrefix()
1415e7058c64SPeter Brune @*/
1416f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[])
1417e7058c64SPeter Brune {
1418e7058c64SPeter Brune   PetscErrorCode ierr;
1419e7058c64SPeter Brune 
1420e7058c64SPeter Brune   PetscFunctionBegin;
1421f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1422e7058c64SPeter Brune   ierr = PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1423e7058c64SPeter Brune   PetscFunctionReturn(0);
1424e7058c64SPeter Brune }
1425e7058c64SPeter Brune 
1426e7058c64SPeter Brune #undef __FUNCT__
1427f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetOptionsPrefix"
1428e7058c64SPeter Brune /*@C
1429f1c6b773SPeter Brune    SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all
1430f1c6b773SPeter Brune    SNESLineSearch options in the database.
1431e7058c64SPeter Brune 
1432e7058c64SPeter Brune    Not Collective
1433e7058c64SPeter Brune 
1434e7058c64SPeter Brune    Input Parameter:
1435cd7522eaSPeter Brune .  linesearch - the SNESLineSearch context
1436e7058c64SPeter Brune 
1437e7058c64SPeter Brune    Output Parameter:
1438e7058c64SPeter Brune .  prefix - pointer to the prefix string used
1439e7058c64SPeter Brune 
14408e557f58SPeter Brune    Notes:
14418e557f58SPeter Brune    On the fortran side, the user should pass in a string 'prefix' of
1442e7058c64SPeter Brune    sufficient length to hold the prefix.
1443e7058c64SPeter Brune 
1444e7058c64SPeter Brune    Level: advanced
1445e7058c64SPeter Brune 
1446f1c6b773SPeter Brune .keywords: SNESLineSearch, get, options, prefix, database
1447e7058c64SPeter Brune 
1448e7058c64SPeter Brune .seealso: SNESAppendOptionsPrefix()
1449e7058c64SPeter Brune @*/
1450f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[])
1451e7058c64SPeter Brune {
1452e7058c64SPeter Brune   PetscErrorCode ierr;
1453e7058c64SPeter Brune 
1454e7058c64SPeter Brune   PetscFunctionBegin;
1455f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1456e7058c64SPeter Brune   ierr = PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1457e7058c64SPeter Brune   PetscFunctionReturn(0);
1458e7058c64SPeter Brune }
1459bf7f4e0aSPeter Brune 
1460bf7f4e0aSPeter Brune #undef __FUNCT__
1461f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetWork"
1462f40b411bSPeter Brune /*@
1463f1c6b773SPeter Brune    SNESLineSearchGetWork - Gets work vectors for the line search.
1464f40b411bSPeter Brune 
1465f40b411bSPeter Brune    Input Parameter:
1466f1c6b773SPeter Brune +  linesearch - the SNESLineSearch context
1467f40b411bSPeter Brune -  nwork - the number of work vectors
1468f40b411bSPeter Brune 
1469f40b411bSPeter Brune    Level: developer
1470f40b411bSPeter Brune 
1471cd7522eaSPeter Brune    Notes:
1472cd7522eaSPeter Brune    This is typically called at the beginning of a SNESLineSearch or SNESLineSearchShell implementation.
1473cd7522eaSPeter Brune 
1474f1c6b773SPeter Brune .keywords: SNESLineSearch, work, vector
1475f40b411bSPeter Brune 
1476f40b411bSPeter Brune .seealso: SNESDefaultGetWork()
1477f40b411bSPeter Brune @*/
1478f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetWork(SNESLineSearch linesearch, PetscInt nwork)
1479bf7f4e0aSPeter Brune {
1480bf7f4e0aSPeter Brune   PetscErrorCode ierr;
1481bf388a1fSBarry Smith 
1482bf7f4e0aSPeter Brune   PetscFunctionBegin;
1483bf7f4e0aSPeter Brune   if (linesearch->vec_sol) {
1484bf7f4e0aSPeter Brune     ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr);
14850ad7597dSKarl Rupp   }
1486*ce94432eSBarry Smith   SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
1487bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1488bf7f4e0aSPeter Brune }
1489bf7f4e0aSPeter Brune 
1490bf7f4e0aSPeter Brune #undef __FUNCT__
1491f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetSuccess"
1492f40b411bSPeter Brune /*@
1493f1c6b773SPeter Brune    SNESLineSearchGetSuccess - Gets the success/failure status of the last line search application
1494f40b411bSPeter Brune 
1495f40b411bSPeter Brune    Input Parameters:
149678bcb3b5SPeter Brune .  linesearch - linesearch context
1497f40b411bSPeter Brune 
1498f40b411bSPeter Brune    Output Parameters:
14998e557f58SPeter Brune .  success - The success or failure status
1500f40b411bSPeter Brune 
1501cd7522eaSPeter Brune    Notes:
1502cd7522eaSPeter Brune    This is typically called after SNESLineSearchApply in order to determine if the line-search failed
1503cd7522eaSPeter Brune    (and set the SNES convergence accordingly).
1504cd7522eaSPeter Brune 
1505f40b411bSPeter Brune    Level: intermediate
1506f40b411bSPeter Brune 
1507f1c6b773SPeter Brune .seealso: SNESLineSearchSetSuccess()
1508f40b411bSPeter Brune @*/
1509f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchGetSuccess(SNESLineSearch linesearch, PetscBool *success)
1510bf7f4e0aSPeter Brune {
1511bf7f4e0aSPeter Brune   PetscFunctionBegin;
1512f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
15136a388c36SPeter Brune   PetscValidPointer(success, 2);
1514f5af7f23SKarl Rupp   if (success) *success = linesearch->success;
1515bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1516bf7f4e0aSPeter Brune }
1517bf7f4e0aSPeter Brune 
1518bf7f4e0aSPeter Brune #undef __FUNCT__
1519f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetSuccess"
1520f40b411bSPeter Brune /*@
1521f1c6b773SPeter Brune    SNESLineSearchSetSuccess - Sets the success/failure status of the last line search application
1522f40b411bSPeter Brune 
1523f40b411bSPeter Brune    Input Parameters:
152478bcb3b5SPeter Brune +  linesearch - linesearch context
15258e557f58SPeter Brune -  success - The success or failure status
1526f40b411bSPeter Brune 
1527cd7522eaSPeter Brune    Notes:
1528cd7522eaSPeter Brune    This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set
1529cd7522eaSPeter Brune    the success or failure of the line search method.
1530cd7522eaSPeter Brune 
153178bcb3b5SPeter Brune    Level: developer
1532f40b411bSPeter Brune 
1533f1c6b773SPeter Brune .seealso: SNESLineSearchGetSuccess()
1534f40b411bSPeter Brune @*/
1535f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchSetSuccess(SNESLineSearch linesearch, PetscBool success)
15366a388c36SPeter Brune {
15376a388c36SPeter Brune   PetscFunctionBegin;
15385fd66863SKarl Rupp   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
15396a388c36SPeter Brune   linesearch->success = success;
15406a388c36SPeter Brune   PetscFunctionReturn(0);
15416a388c36SPeter Brune }
15426a388c36SPeter Brune 
15436a388c36SPeter Brune #undef __FUNCT__
1544f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetVIFunctions"
15459bd66eb0SPeter Brune /*@C
1546f1c6b773SPeter Brune    SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.
15479bd66eb0SPeter Brune 
15489bd66eb0SPeter Brune    Input Parameters:
15499bd66eb0SPeter Brune +  snes - nonlinear context obtained from SNESCreate()
15509bd66eb0SPeter Brune .  projectfunc - function for projecting the function to the bounds
15519bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
15529bd66eb0SPeter Brune 
15539bd66eb0SPeter Brune    Logically Collective on SNES
15549bd66eb0SPeter Brune 
15559bd66eb0SPeter Brune    Calling sequence of projectfunc:
15569bd66eb0SPeter Brune .vb
15579bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X)
15589bd66eb0SPeter Brune .ve
15599bd66eb0SPeter Brune 
15609bd66eb0SPeter Brune     Input parameters for projectfunc:
15619bd66eb0SPeter Brune +   snes - nonlinear context
15629bd66eb0SPeter Brune -   X - current solution
15639bd66eb0SPeter Brune 
1564cd7522eaSPeter Brune     Output parameters for projectfunc:
15659bd66eb0SPeter Brune .   X - Projected solution
15669bd66eb0SPeter Brune 
15679bd66eb0SPeter Brune    Calling sequence of normfunc:
15689bd66eb0SPeter Brune .vb
15699bd66eb0SPeter Brune    projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm)
15709bd66eb0SPeter Brune .ve
15719bd66eb0SPeter Brune 
1572cd7522eaSPeter Brune     Input parameters for normfunc:
15739bd66eb0SPeter Brune +   snes - nonlinear context
15749bd66eb0SPeter Brune .   X - current solution
15759bd66eb0SPeter Brune -   F - current residual
15769bd66eb0SPeter Brune 
1577cd7522eaSPeter Brune     Output parameters for normfunc:
15789bd66eb0SPeter Brune .   fnorm - VI-specific norm of the function
15799bd66eb0SPeter Brune 
1580cd7522eaSPeter Brune     Notes:
1581cd7522eaSPeter Brune     The VI solvers require projection of the solution to the feasible set.  projectfunc should implement this.
1582cd7522eaSPeter Brune 
1583cd7522eaSPeter Brune     The VI solvers require special evaluation of the function norm such that the norm is only calculated
1584cd7522eaSPeter Brune     on the inactive set.  This should be implemented by normfunc.
15859bd66eb0SPeter Brune 
15869bd66eb0SPeter Brune     Level: developer
15879bd66eb0SPeter Brune 
15889bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, set, line search
15899bd66eb0SPeter Brune 
1590f1c6b773SPeter Brune .seealso: SNESLineSearchGetVIFunctions(), SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck()
15919bd66eb0SPeter Brune @*/
1592f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc)
15939bd66eb0SPeter Brune {
15949bd66eb0SPeter Brune   PetscFunctionBegin;
1595f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
15969bd66eb0SPeter Brune   if (projectfunc) linesearch->ops->viproject = projectfunc;
15979bd66eb0SPeter Brune   if (normfunc) linesearch->ops->vinorm = normfunc;
15989bd66eb0SPeter Brune   PetscFunctionReturn(0);
15999bd66eb0SPeter Brune }
16009bd66eb0SPeter Brune 
16019bd66eb0SPeter Brune #undef __FUNCT__
1602f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchGetVIFunctions"
16039bd66eb0SPeter Brune /*@C
1604f1c6b773SPeter Brune    SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.
16059bd66eb0SPeter Brune 
16069bd66eb0SPeter Brune    Input Parameters:
16079bd66eb0SPeter Brune .  snes - nonlinear context obtained from SNESCreate()
16089bd66eb0SPeter Brune 
16099bd66eb0SPeter Brune    Output Parameters:
16109bd66eb0SPeter Brune +  projectfunc - function for projecting the function to the bounds
16119bd66eb0SPeter Brune -  normfunc - function for computing the norm of an active set
16129bd66eb0SPeter Brune 
16139bd66eb0SPeter Brune    Logically Collective on SNES
16149bd66eb0SPeter Brune 
16159bd66eb0SPeter Brune     Level: developer
16169bd66eb0SPeter Brune 
16179bd66eb0SPeter Brune .keywords: SNES, line search, VI, nonlinear, get, line search
16189bd66eb0SPeter Brune 
1619f1c6b773SPeter Brune .seealso: SNESLineSearchSetVIFunctions(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
16209bd66eb0SPeter Brune @*/
1621f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc)
16229bd66eb0SPeter Brune {
16239bd66eb0SPeter Brune   PetscFunctionBegin;
16249bd66eb0SPeter Brune   if (projectfunc) *projectfunc = linesearch->ops->viproject;
16259bd66eb0SPeter Brune   if (normfunc) *normfunc = linesearch->ops->vinorm;
16269bd66eb0SPeter Brune   PetscFunctionReturn(0);
16279bd66eb0SPeter Brune }
16289bd66eb0SPeter Brune 
16299bd66eb0SPeter Brune #undef __FUNCT__
1630f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchRegister"
1631bf7f4e0aSPeter Brune /*@C
1632f1c6b773SPeter Brune   SNESLineSearchRegister - See SNESLineSearchRegisterDynamic()
1633bf7f4e0aSPeter Brune 
1634bf7f4e0aSPeter Brune   Level: advanced
1635bf7f4e0aSPeter Brune @*/
1636f1c6b773SPeter Brune PetscErrorCode  SNESLineSearchRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNESLineSearch))
1637bf7f4e0aSPeter Brune {
1638bf7f4e0aSPeter Brune   char           fullname[PETSC_MAX_PATH_LEN];
1639bf7f4e0aSPeter Brune   PetscErrorCode ierr;
1640bf7f4e0aSPeter Brune 
1641bf7f4e0aSPeter Brune   PetscFunctionBegin;
1642140e18c1SBarry Smith   ierr = PetscFunctionListConcat(path,name,fullname);CHKERRQ(ierr);
1643140e18c1SBarry Smith   ierr = PetscFunctionListAdd(PETSC_COMM_WORLD,&SNESLineSearchList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
1644bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
1645bf7f4e0aSPeter Brune }
1646