xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision 8c5add6abd9f0e97420a8d072677ee1e39293c49)
1af0996ceSBarry Smith #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;
757a83faaSBarry Smith PetscLogEvent SNESLINESEARCH_Apply;
8bf7f4e0aSPeter Brune 
9dcf2fd19SBarry Smith /*@
10f6dfbefdSBarry Smith   SNESLineSearchMonitorCancel - Clears all the monitor functions for a `SNESLineSearch` object.
11dcf2fd19SBarry Smith 
12c3339decSBarry Smith   Logically Collective
13dcf2fd19SBarry Smith 
142fe279fdSBarry Smith   Input Parameter:
15f6dfbefdSBarry Smith . ls - the `SNESLineSearch` context
16dcf2fd19SBarry Smith 
17dcf2fd19SBarry Smith   Options Database Key:
18dcf2fd19SBarry Smith . -snes_linesearch_monitor_cancel - cancels all monitors that have been hardwired
19f6dfbefdSBarry Smith     into a code by calls to `SNESLineSearchMonitorSet()`, but does not cancel those
20dcf2fd19SBarry Smith     set via the options database
21dcf2fd19SBarry Smith 
2220f4b53cSBarry Smith   Level: advanced
2320f4b53cSBarry Smith 
24dcf2fd19SBarry Smith   Notes:
25f6dfbefdSBarry Smith   There is no way to clear one specific monitor from a `SNESLineSearch` object.
26dcf2fd19SBarry Smith 
27420bcc1bSBarry Smith   This does not clear the monitor set with `SNESLineSearchSetDefaultMonitor()` use `SNESLineSearchSetDefaultMonitor`(`ls`,`NULL`) to cancel it
28dcf2fd19SBarry Smith   that one.
29dcf2fd19SBarry Smith 
30420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorDefault()`, `SNESLineSearchMonitorSet()`
31dcf2fd19SBarry Smith @*/
32d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchMonitorCancel(SNESLineSearch ls)
33d71ae5a4SJacob Faibussowitsch {
34dcf2fd19SBarry Smith   PetscInt i;
35dcf2fd19SBarry Smith 
36dcf2fd19SBarry Smith   PetscFunctionBegin;
37dcf2fd19SBarry Smith   PetscValidHeaderSpecific(ls, SNESLINESEARCH_CLASSID, 1);
38dcf2fd19SBarry Smith   for (i = 0; i < ls->numbermonitors; i++) {
3948a46eb9SPierre Jolivet     if (ls->monitordestroy[i]) PetscCall((*ls->monitordestroy[i])(&ls->monitorcontext[i]));
40dcf2fd19SBarry Smith   }
41dcf2fd19SBarry Smith   ls->numbermonitors = 0;
423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
43dcf2fd19SBarry Smith }
44dcf2fd19SBarry Smith 
45dcf2fd19SBarry Smith /*@
46dcf2fd19SBarry Smith   SNESLineSearchMonitor - runs the user provided monitor routines, if they exist
47dcf2fd19SBarry Smith 
48c3339decSBarry Smith   Collective
49dcf2fd19SBarry Smith 
502fe279fdSBarry Smith   Input Parameter:
51dcf2fd19SBarry Smith . ls - the linesearch object
52dcf2fd19SBarry Smith 
532fe279fdSBarry Smith   Level: developer
542fe279fdSBarry Smith 
55f6dfbefdSBarry Smith   Note:
56420bcc1bSBarry Smith   This routine is called by the `SNESLineSearch` implementations.
57dcf2fd19SBarry Smith   It does not typically need to be called by the user.
58dcf2fd19SBarry Smith 
59420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorSet()`
60dcf2fd19SBarry Smith @*/
61d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchMonitor(SNESLineSearch ls)
62d71ae5a4SJacob Faibussowitsch {
63dcf2fd19SBarry Smith   PetscInt i, n = ls->numbermonitors;
64dcf2fd19SBarry Smith 
65dcf2fd19SBarry Smith   PetscFunctionBegin;
6648a46eb9SPierre Jolivet   for (i = 0; i < n; i++) PetscCall((*ls->monitorftns[i])(ls, ls->monitorcontext[i]));
673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
68dcf2fd19SBarry Smith }
69dcf2fd19SBarry Smith 
70dcf2fd19SBarry Smith /*@C
71dcf2fd19SBarry Smith   SNESLineSearchMonitorSet - Sets an ADDITIONAL function that is to be used at every
72dcf2fd19SBarry Smith   iteration of the nonlinear solver to display the iteration's
73dcf2fd19SBarry Smith   progress.
74dcf2fd19SBarry Smith 
75c3339decSBarry Smith   Logically Collective
76dcf2fd19SBarry Smith 
77dcf2fd19SBarry Smith   Input Parameters:
78f6dfbefdSBarry Smith + ls             - the `SNESLineSearch` context
79dcf2fd19SBarry Smith . f              - the monitor function
80420bcc1bSBarry Smith . mctx           - [optional] user-defined context for private data for the monitor routine (use `NULL` if no context is desired)
8149abdd8aSBarry Smith - monitordestroy - [optional] routine that frees monitor context (may be `NULL`), see `PetscCtxDestroyFn` for the calling sequence
82420bcc1bSBarry Smith 
83420bcc1bSBarry Smith   Calling sequence of `f`:
84420bcc1bSBarry Smith + ls   - the `SNESLineSearch` context
85420bcc1bSBarry Smith - mctx - [optional] user-defined context for private data for the monitor routine
86420bcc1bSBarry Smith 
8720f4b53cSBarry Smith   Level: intermediate
88dcf2fd19SBarry Smith 
89f6dfbefdSBarry Smith   Note:
90dcf2fd19SBarry Smith   Several different monitoring routines may be set by calling
91f6dfbefdSBarry Smith   `SNESLineSearchMonitorSet()` multiple times; all will be called in the
92dcf2fd19SBarry Smith   order in which they were set.
93dcf2fd19SBarry Smith 
94420bcc1bSBarry Smith   Fortran Note:
95f6dfbefdSBarry Smith   Only a single monitor function can be set for each `SNESLineSearch` object
96dcf2fd19SBarry Smith 
9749abdd8aSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorDefault()`, `SNESLineSearchMonitorCancel()`, `PetscCtxDestroyFn`
98dcf2fd19SBarry Smith @*/
9949abdd8aSBarry Smith PetscErrorCode SNESLineSearchMonitorSet(SNESLineSearch ls, PetscErrorCode (*f)(SNESLineSearch ls, void *mctx), void *mctx, PetscCtxDestroyFn *monitordestroy)
100d71ae5a4SJacob Faibussowitsch {
101dcf2fd19SBarry Smith   PetscFunctionBegin;
102dcf2fd19SBarry Smith   PetscValidHeaderSpecific(ls, SNESLINESEARCH_CLASSID, 1);
103453a69bbSBarry Smith   for (PetscInt i = 0; i < ls->numbermonitors; i++) {
104453a69bbSBarry Smith     PetscBool identical;
105453a69bbSBarry Smith 
106453a69bbSBarry Smith     PetscCall(PetscMonitorCompare((PetscErrorCode (*)(void))(PetscVoidFn *)f, mctx, monitordestroy, (PetscErrorCode (*)(void))(PetscVoidFn *)ls->monitorftns[i], ls->monitorcontext[i], ls->monitordestroy[i], &identical));
1073ba16761SJacob Faibussowitsch     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
10878064530SBarry Smith   }
1095f80ce2aSJacob Faibussowitsch   PetscCheck(ls->numbermonitors < MAXSNESLSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many monitors set");
110dcf2fd19SBarry Smith   ls->monitorftns[ls->numbermonitors]      = f;
111dcf2fd19SBarry Smith   ls->monitordestroy[ls->numbermonitors]   = monitordestroy;
112835f2295SStefano Zampini   ls->monitorcontext[ls->numbermonitors++] = mctx;
1133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
114dcf2fd19SBarry Smith }
115dcf2fd19SBarry Smith 
116dcf2fd19SBarry Smith /*@C
117f6dfbefdSBarry Smith   SNESLineSearchMonitorSolutionUpdate - Monitors each update of the function value the linesearch tries
118dcf2fd19SBarry Smith 
119c3339decSBarry Smith   Collective
120dcf2fd19SBarry Smith 
121dcf2fd19SBarry Smith   Input Parameters:
122420bcc1bSBarry Smith + ls - the `SNESLineSearch` object
123f6dfbefdSBarry Smith - vf - the context for the monitor, in this case it is an `PetscViewerAndFormat`
124dcf2fd19SBarry Smith 
125420bcc1bSBarry Smith   Options Database Key:
126420bcc1bSBarry Smith . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine
127dcf2fd19SBarry Smith 
128420bcc1bSBarry Smith   Level: developer
129420bcc1bSBarry Smith 
130420bcc1bSBarry Smith   This is not normally called directly but is passed to `SNESLineSearchMonitorSet()`
131420bcc1bSBarry Smith 
132420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorSet()`, `SNESMonitorSolution()`
133dcf2fd19SBarry Smith @*/
134d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchMonitorSolutionUpdate(SNESLineSearch ls, PetscViewerAndFormat *vf)
135d71ae5a4SJacob Faibussowitsch {
136d12e167eSBarry Smith   PetscViewer viewer = vf->viewer;
137dcf2fd19SBarry Smith   Vec         Y, W, G;
138dcf2fd19SBarry Smith 
139dcf2fd19SBarry Smith   PetscFunctionBegin;
1409566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetVecs(ls, NULL, NULL, &Y, &W, &G));
1419566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(viewer, vf->format));
1429566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted update to solution \n"));
1439566063dSJacob Faibussowitsch   PetscCall(VecView(Y, viewer));
1449566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted new solution \n"));
1459566063dSJacob Faibussowitsch   PetscCall(VecView(W, viewer));
1469566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted updated function value\n"));
1479566063dSJacob Faibussowitsch   PetscCall(VecView(G, viewer));
1489566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(viewer));
1493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
150dcf2fd19SBarry Smith }
151dcf2fd19SBarry Smith 
152f40b411bSPeter Brune /*@
153420bcc1bSBarry Smith   SNESLineSearchCreate - Creates a `SNESLineSearch` context.
154f40b411bSPeter Brune 
155f6dfbefdSBarry Smith   Logically Collective
156f40b411bSPeter Brune 
1572fe279fdSBarry Smith   Input Parameter:
158f6dfbefdSBarry Smith . comm - MPI communicator for the line search (typically from the associated `SNES` context).
159f40b411bSPeter Brune 
1602fe279fdSBarry Smith   Output Parameter:
1618e557f58SPeter Brune . outlinesearch - the new line search context
162f40b411bSPeter Brune 
163162e0bf5SPeter Brune   Level: developer
164162e0bf5SPeter Brune 
165f6dfbefdSBarry Smith   Note:
166420bcc1bSBarry Smith   The preferred calling sequence is to use `SNESGetLineSearch()` to acquire the `SNESLineSearch` instance
167f6dfbefdSBarry Smith   already associated with the `SNES`.
168f40b411bSPeter Brune 
169420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `LineSearchDestroy()`, `SNESGetLineSearch()`
170f40b411bSPeter Brune @*/
171d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch)
172d71ae5a4SJacob Faibussowitsch {
173f1c6b773SPeter Brune   SNESLineSearch linesearch;
174bf388a1fSBarry Smith 
175bf7f4e0aSPeter Brune   PetscFunctionBegin;
1764f572ea9SToby Isaac   PetscAssertPointer(outlinesearch, 2);
1779566063dSJacob Faibussowitsch   PetscCall(SNESInitializePackage());
178f5af7f23SKarl Rupp 
1799566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(linesearch, SNESLINESEARCH_CLASSID, "SNESLineSearch", "Linesearch", "SNESLineSearch", comm, SNESLineSearchDestroy, SNESLineSearchView));
1800298fd71SBarry Smith   linesearch->vec_sol_new  = NULL;
1810298fd71SBarry Smith   linesearch->vec_func_new = NULL;
1820298fd71SBarry Smith   linesearch->vec_sol      = NULL;
1830298fd71SBarry Smith   linesearch->vec_func     = NULL;
1840298fd71SBarry Smith   linesearch->vec_update   = NULL;
1859bd66eb0SPeter Brune 
186bf7f4e0aSPeter Brune   linesearch->lambda       = 1.0;
187bf7f4e0aSPeter Brune   linesearch->fnorm        = 1.0;
188bf7f4e0aSPeter Brune   linesearch->ynorm        = 1.0;
189bf7f4e0aSPeter Brune   linesearch->xnorm        = 1.0;
190422a814eSBarry Smith   linesearch->result       = SNES_LINESEARCH_SUCCEEDED;
191bf7f4e0aSPeter Brune   linesearch->norms        = PETSC_TRUE;
192bf7f4e0aSPeter Brune   linesearch->keeplambda   = PETSC_FALSE;
193bf7f4e0aSPeter Brune   linesearch->damping      = 1.0;
194a99ef635SJonas Heinzmann   linesearch->maxlambda    = 1.0;
195a99ef635SJonas Heinzmann   linesearch->minlambda    = 1e-12;
196516fe3c3SPeter Brune   linesearch->rtol         = 1e-8;
197516fe3c3SPeter Brune   linesearch->atol         = 1e-15;
198516fe3c3SPeter Brune   linesearch->ltol         = 1e-8;
1990298fd71SBarry Smith   linesearch->precheckctx  = NULL;
2000298fd71SBarry Smith   linesearch->postcheckctx = NULL;
201a99ef635SJonas Heinzmann   linesearch->max_it       = 1;
202bf7f4e0aSPeter Brune   linesearch->setupcalled  = PETSC_FALSE;
2033add74b1SBarry Smith   linesearch->monitor      = NULL;
204bf7f4e0aSPeter Brune   *outlinesearch           = linesearch;
2053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
206bf7f4e0aSPeter Brune }
207bf7f4e0aSPeter Brune 
208f40b411bSPeter Brune /*@
20978bcb3b5SPeter Brune   SNESLineSearchSetUp - Prepares the line search for being applied by allocating
21078bcb3b5SPeter Brune   any required vectors.
211f40b411bSPeter Brune 
212c3339decSBarry Smith   Collective
213f40b411bSPeter Brune 
2142fe279fdSBarry Smith   Input Parameter:
215f6dfbefdSBarry Smith . linesearch - The `SNESLineSearch` instance.
216f40b411bSPeter Brune 
21720f4b53cSBarry Smith   Level: advanced
21820f4b53cSBarry Smith 
219f6dfbefdSBarry Smith   Note:
220f6dfbefdSBarry Smith   For most cases, this needn't be called by users or outside of `SNESLineSearchApply()`.
221cd7522eaSPeter Brune   The only current case where this is called outside of this is for the VI
22278bcb3b5SPeter Brune   solvers, which modify the solution and work vectors before the first call
223f6dfbefdSBarry Smith   of `SNESLineSearchApply()`, requiring the `SNESLineSearch` work vectors to be
224cd7522eaSPeter Brune   allocated upfront.
225cd7522eaSPeter Brune 
226420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchReset()`
227f40b411bSPeter Brune @*/
228d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch)
229d71ae5a4SJacob Faibussowitsch {
230bf7f4e0aSPeter Brune   PetscFunctionBegin;
23148a46eb9SPierre Jolivet   if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC));
232bf7f4e0aSPeter Brune   if (!linesearch->setupcalled) {
23348a46eb9SPierre Jolivet     if (!linesearch->vec_sol_new) PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new));
23448a46eb9SPierre Jolivet     if (!linesearch->vec_func_new) PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new));
235213acdd3SPierre Jolivet     PetscTryTypeMethod(linesearch, setup);
2369566063dSJacob Faibussowitsch     if (!linesearch->ops->snesfunc) PetscCall(SNESLineSearchSetFunction(linesearch, SNESComputeFunction));
237bf7f4e0aSPeter Brune     linesearch->lambda      = linesearch->damping;
238bf7f4e0aSPeter Brune     linesearch->setupcalled = PETSC_TRUE;
239bf7f4e0aSPeter Brune   }
2403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
241bf7f4e0aSPeter Brune }
242bf7f4e0aSPeter Brune 
243f40b411bSPeter Brune /*@
244f6dfbefdSBarry Smith   SNESLineSearchReset - Undoes the `SNESLineSearchSetUp()` and deletes any `Vec`s or `Mat`s allocated by the line search.
245f40b411bSPeter Brune 
246c3339decSBarry Smith   Collective
247f40b411bSPeter Brune 
2482fe279fdSBarry Smith   Input Parameter:
249f6dfbefdSBarry Smith . linesearch - The `SNESLineSearch` instance.
250f40b411bSPeter Brune 
25120f4b53cSBarry Smith   Level: developer
25220f4b53cSBarry Smith 
253f6dfbefdSBarry Smith   Note:
254f6dfbefdSBarry Smith   Usually only called by `SNESReset()`
255f190f2fcSBarry Smith 
256420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchSetUp()`
257f40b411bSPeter Brune @*/
258d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch)
259d71ae5a4SJacob Faibussowitsch {
260bf7f4e0aSPeter Brune   PetscFunctionBegin;
261213acdd3SPierre Jolivet   PetscTryTypeMethod(linesearch, reset);
262f5af7f23SKarl Rupp 
2639566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&linesearch->vec_sol_new));
2649566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&linesearch->vec_func_new));
265bf7f4e0aSPeter Brune 
2669566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(linesearch->nwork, &linesearch->work));
267f5af7f23SKarl Rupp 
268bf7f4e0aSPeter Brune   linesearch->nwork       = 0;
269bf7f4e0aSPeter Brune   linesearch->setupcalled = PETSC_FALSE;
2703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
271bf7f4e0aSPeter Brune }
272bf7f4e0aSPeter Brune 
273ed07d7d7SPeter Brune /*@C
274f6dfbefdSBarry Smith   SNESLineSearchSetFunction - Sets the function evaluation used by the `SNES` line search
275f6dfbefdSBarry Smith   `
276e4094ef1SJacob Faibussowitsch 
277ed07d7d7SPeter Brune   Input Parameters:
278e4094ef1SJacob Faibussowitsch + linesearch - the `SNESLineSearch` context
279e4094ef1SJacob Faibussowitsch - func       - function evaluation routine, this is usually the function provided with `SNESSetFunction()`
280ed07d7d7SPeter Brune 
281420bcc1bSBarry Smith   Calling sequence of `func`:
282420bcc1bSBarry Smith + snes - the `SNES` with which the `SNESLineSearch` context is associated with
283420bcc1bSBarry Smith . x    - the input vector
284420bcc1bSBarry Smith - f    - the computed value of the function
285420bcc1bSBarry Smith 
286ed07d7d7SPeter Brune   Level: developer
287ed07d7d7SPeter Brune 
288420bcc1bSBarry Smith   Note:
289420bcc1bSBarry Smith   By default the `SNESLineSearch` uses the function provided by `SNESSetFunction()` so this is rarely needed
290420bcc1bSBarry Smith 
291420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESSetFunction()`
292ed07d7d7SPeter Brune @*/
293420bcc1bSBarry Smith PetscErrorCode SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES snes, Vec x, Vec f))
294d71ae5a4SJacob Faibussowitsch {
295ed07d7d7SPeter Brune   PetscFunctionBegin;
296ed07d7d7SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
297ed07d7d7SPeter Brune   linesearch->ops->snesfunc = func;
2983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
299ed07d7d7SPeter Brune }
300ed07d7d7SPeter Brune 
30186d74e61SPeter Brune /*@C
302420bcc1bSBarry Smith   SNESLineSearchSetPreCheck - Sets a function that is called after the initial search direction has been computed but
303420bcc1bSBarry Smith   before the line search routine has been applied. Allows adjusting the result of (usually a linear solve) that
304f190f2fcSBarry Smith   determined the search direction.
30586d74e61SPeter Brune 
306c3339decSBarry Smith   Logically Collective
30786d74e61SPeter Brune 
30886d74e61SPeter Brune   Input Parameters:
309f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context
310420bcc1bSBarry Smith . func       - [optional] function evaluation routine
31120f4b53cSBarry Smith - ctx        - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
31286d74e61SPeter Brune 
313420bcc1bSBarry Smith   Calling sequence of `func`:
314420bcc1bSBarry Smith + ls        - the `SNESLineSearch` context
315420bcc1bSBarry Smith . x         - the current solution
316a678f235SPierre Jolivet . d         - the current search direction
317420bcc1bSBarry Smith . changed_d - indicates if the search direction has been changed
318420bcc1bSBarry Smith - ctx       - the context passed to `SNESLineSearchSetPreCheck()`
319420bcc1bSBarry Smith 
32086d74e61SPeter Brune   Level: intermediate
32186d74e61SPeter Brune 
322f6dfbefdSBarry Smith   Note:
323420bcc1bSBarry Smith   Use `SNESLineSearchSetPostCheck()` to change the step after the line search is complete.
324f0b84518SBarry Smith 
325f0b84518SBarry Smith   Use `SNESVISetVariableBounds()` and `SNESVISetComputeVariableBounds()` to cause `SNES` to automatically control the ranges of variables allowed.
326f0b84518SBarry Smith 
327420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()`,
328869ba2dcSBarry Smith           `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`
329f0b84518SBarry Smith 
33086d74e61SPeter Brune @*/
331420bcc1bSBarry Smith PetscErrorCode SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch ls, Vec x, Vec d, PetscBool *changed_d, void *ctx), void *ctx)
332d71ae5a4SJacob Faibussowitsch {
3339bd66eb0SPeter Brune   PetscFunctionBegin;
334f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
335f190f2fcSBarry Smith   if (func) linesearch->ops->precheck = func;
33686d74e61SPeter Brune   if (ctx) linesearch->precheckctx = ctx;
3373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
33886d74e61SPeter Brune }
33986d74e61SPeter Brune 
34086d74e61SPeter Brune /*@C
34153deb899SBarry Smith   SNESLineSearchGetPreCheck - Gets the pre-check function for the line search routine.
34286d74e61SPeter Brune 
343f899ff85SJose E. Roman   Input Parameter:
344f6dfbefdSBarry Smith . linesearch - the `SNESLineSearch` context
34586d74e61SPeter Brune 
34686d74e61SPeter Brune   Output Parameters:
347420bcc1bSBarry Smith + func - [optional] function evaluation routine,  for calling sequence see `SNESLineSearchSetPreCheck()`
34820f4b53cSBarry Smith - ctx  - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
34986d74e61SPeter Brune 
35086d74e61SPeter Brune   Level: intermediate
35186d74e61SPeter Brune 
352420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchSetPostCheck()`
35386d74e61SPeter Brune @*/
354d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch, Vec, Vec, PetscBool *, void *), void **ctx)
355d71ae5a4SJacob Faibussowitsch {
3569bd66eb0SPeter Brune   PetscFunctionBegin;
357f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
358f190f2fcSBarry Smith   if (func) *func = linesearch->ops->precheck;
35986d74e61SPeter Brune   if (ctx) *ctx = linesearch->precheckctx;
3603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
36186d74e61SPeter Brune }
36286d74e61SPeter Brune 
36386d74e61SPeter Brune /*@C
364f190f2fcSBarry Smith   SNESLineSearchSetPostCheck - Sets a user function that is called after the line search has been applied to determine the step
365f190f2fcSBarry Smith   direction and length. Allows the user a chance to change or override the decision of the line search routine
36686d74e61SPeter Brune 
36720f4b53cSBarry Smith   Logically Collective
36886d74e61SPeter Brune 
36986d74e61SPeter Brune   Input Parameters:
370f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context
371420bcc1bSBarry Smith . func       - [optional] function evaluation routine
37220f4b53cSBarry Smith - ctx        - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
37386d74e61SPeter Brune 
374420bcc1bSBarry Smith   Calling sequence of `func`:
375420bcc1bSBarry Smith + ls        - the `SNESLineSearch` context
376420bcc1bSBarry Smith . x         - the current solution
377a678f235SPierre Jolivet . d         - the current search direction
378a678f235SPierre Jolivet . w         - $ w = x + lambda*d $ for some lambda
379420bcc1bSBarry Smith . changed_d - indicates if the search direction `d` has been changed
380420bcc1bSBarry Smith . changed_w - indicates `w` has been changed
381420bcc1bSBarry Smith - ctx       - the context passed to `SNESLineSearchSetPreCheck()`
382420bcc1bSBarry Smith 
38386d74e61SPeter Brune   Level: intermediate
38486d74e61SPeter Brune 
385f0b84518SBarry Smith   Notes:
3866b095a85SStefano Zampini   Use `SNESLineSearchSetPreCheck()` to change the step before the line search is completed.
3876b095a85SStefano Zampini   The calling sequence of the callback does not contain the current scaling factor. To access the value, use `SNESLineSearchGetLambda()`.
388f0b84518SBarry Smith 
389f0b84518SBarry Smith   Use `SNESVISetVariableBounds()` and `SNESVISetComputeVariableBounds()` to cause `SNES` to automatically control the ranges of variables allowed.
390f0b84518SBarry Smith 
391420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchGetPostCheck()`,
392e4094ef1SJacob Faibussowitsch           `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`
39386d74e61SPeter Brune @*/
394420bcc1bSBarry Smith PetscErrorCode SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch ls, Vec x, Vec d, Vec w, PetscBool *changed_d, PetscBool *changed_w, void *ctx), void *ctx)
395d71ae5a4SJacob Faibussowitsch {
39686d74e61SPeter Brune   PetscFunctionBegin;
397f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
398f190f2fcSBarry Smith   if (func) linesearch->ops->postcheck = func;
39986d74e61SPeter Brune   if (ctx) linesearch->postcheckctx = ctx;
4003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
40186d74e61SPeter Brune }
40286d74e61SPeter Brune 
40386d74e61SPeter Brune /*@C
404f1c6b773SPeter Brune   SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine.
40586d74e61SPeter Brune 
406f899ff85SJose E. Roman   Input Parameter:
407f6dfbefdSBarry Smith . linesearch - the `SNESLineSearch` context
40886d74e61SPeter Brune 
40986d74e61SPeter Brune   Output Parameters:
410420bcc1bSBarry Smith + func - [optional] function evaluation routine, see for the calling sequence `SNESLineSearchSetPostCheck()`
41120f4b53cSBarry Smith - ctx  - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
41286d74e61SPeter Brune 
41386d74e61SPeter Brune   Level: intermediate
41486d74e61SPeter Brune 
415420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`
41686d74e61SPeter Brune @*/
417d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx)
418d71ae5a4SJacob Faibussowitsch {
4199bd66eb0SPeter Brune   PetscFunctionBegin;
420f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
421f190f2fcSBarry Smith   if (func) *func = linesearch->ops->postcheck;
42286d74e61SPeter Brune   if (ctx) *ctx = linesearch->postcheckctx;
4233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
42486d74e61SPeter Brune }
42586d74e61SPeter Brune 
426f40b411bSPeter Brune /*@
427f1c6b773SPeter Brune   SNESLineSearchPreCheck - Prepares the line search for being applied.
428f40b411bSPeter Brune 
429c3339decSBarry Smith   Logically Collective
430f40b411bSPeter Brune 
431f40b411bSPeter Brune   Input Parameters:
4327b1df9c1SPeter Brune + linesearch - The linesearch instance.
4337b1df9c1SPeter Brune . X          - The current solution
4347b1df9c1SPeter Brune - Y          - The step direction
435f40b411bSPeter Brune 
4362fe279fdSBarry Smith   Output Parameter:
437420bcc1bSBarry Smith . changed - Indicator that the precheck routine has changed `Y`
438f40b411bSPeter Brune 
439f0b84518SBarry Smith   Level: advanced
440f40b411bSPeter Brune 
441f6dfbefdSBarry Smith   Note:
442420bcc1bSBarry Smith   This calls any function provided with `SNESLineSearchSetPreCheck()` and is called automatically inside the line search routines
443f6dfbefdSBarry Smith 
444420bcc1bSBarry Smith   Developer Note:
445420bcc1bSBarry Smith   The use of `PetscObjectGetState()` would eliminate the need for the `changed` argument to be provided
446420bcc1bSBarry Smith 
447420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`,
448fe8e7dddSPierre Jolivet           `SNESLineSearchGetPostCheck()`
449f40b411bSPeter Brune @*/
450d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed)
451d71ae5a4SJacob Faibussowitsch {
452bf7f4e0aSPeter Brune   PetscFunctionBegin;
453bf7f4e0aSPeter Brune   *changed = PETSC_FALSE;
4546b2b7091SBarry Smith   if (linesearch->ops->precheck) {
455dbbe0bcdSBarry Smith     PetscUseTypeMethod(linesearch, precheck, X, Y, changed, linesearch->precheckctx);
45638bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch, *changed, 4);
457bf7f4e0aSPeter Brune   }
4583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
459bf7f4e0aSPeter Brune }
460bf7f4e0aSPeter Brune 
461f40b411bSPeter Brune /*@
462ef46b1a6SStefano Zampini   SNESLineSearchPostCheck - Hook to modify step direction or updated solution after a successful linesearch
463f40b411bSPeter Brune 
464c3339decSBarry Smith   Logically Collective
465f40b411bSPeter Brune 
466f40b411bSPeter Brune   Input Parameters:
4677b1df9c1SPeter Brune + linesearch - The line search context
4687b1df9c1SPeter Brune . X          - The last solution
4697b1df9c1SPeter Brune . Y          - The step direction
4706b095a85SStefano Zampini - W          - The updated solution, `W = X - lambda * Y` for some lambda
471f40b411bSPeter Brune 
472f40b411bSPeter Brune   Output Parameters:
4736b095a85SStefano Zampini + changed_Y - Indicator if the direction `Y` has been changed.
474420bcc1bSBarry Smith - changed_W - Indicator if the new candidate solution `W` has been changed.
475f40b411bSPeter Brune 
47620f4b53cSBarry Smith   Level: developer
47720f4b53cSBarry Smith 
478f6dfbefdSBarry Smith   Note:
479420bcc1bSBarry Smith   This calls any function provided with `SNESLineSearchSetPostCheck()` and is called automatically inside the line search routines
480f6dfbefdSBarry Smith 
481420bcc1bSBarry Smith   Developer Note:
482420bcc1bSBarry Smith   The use of `PetscObjectGetState()` would eliminate the need for the `changed_Y` and `changed_W` arguments to be provided
483420bcc1bSBarry Smith 
484420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPrecheck()`, `SNESLineSearchGetPrecheck()`
485f40b411bSPeter Brune @*/
486d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W)
487d71ae5a4SJacob Faibussowitsch {
488bf7f4e0aSPeter Brune   PetscFunctionBegin;
489bf7f4e0aSPeter Brune   *changed_Y = PETSC_FALSE;
490bf7f4e0aSPeter Brune   *changed_W = PETSC_FALSE;
4916b2b7091SBarry Smith   if (linesearch->ops->postcheck) {
492dbbe0bcdSBarry Smith     PetscUseTypeMethod(linesearch, postcheck, X, Y, W, changed_Y, changed_W, linesearch->postcheckctx);
49338bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch, *changed_Y, 5);
49438bcdd5aSPeter Brune     PetscValidLogicalCollectiveBool(linesearch, *changed_W, 6);
49586d74e61SPeter Brune   }
4963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
49786d74e61SPeter Brune }
49886d74e61SPeter Brune 
49986d74e61SPeter Brune /*@C
5001d27aa22SBarry Smith   SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration {cite}`hindmarsh1996time`
50186d74e61SPeter Brune 
502c3339decSBarry Smith   Logically Collective
50386d74e61SPeter Brune 
5044165533cSJose E. Roman   Input Parameters:
505420bcc1bSBarry Smith + linesearch - the line search context
50686d74e61SPeter Brune . X          - base state for this step
507907376e6SBarry Smith - ctx        - context for this function
50886d74e61SPeter Brune 
50997bb3fdcSJose E. Roman   Input/Output Parameter:
51097bb3fdcSJose E. Roman . Y - correction, possibly modified
51197bb3fdcSJose E. Roman 
51297bb3fdcSJose E. Roman   Output Parameter:
513420bcc1bSBarry Smith . changed - flag indicating that `Y` was modified
51486d74e61SPeter Brune 
515420bcc1bSBarry Smith   Options Database Keys:
516cd7522eaSPeter Brune + -snes_linesearch_precheck_picard       - activate this routine
517cd7522eaSPeter Brune - -snes_linesearch_precheck_picard_angle - angle
51886d74e61SPeter Brune 
51986d74e61SPeter Brune   Level: advanced
52086d74e61SPeter Brune 
52186d74e61SPeter Brune   Notes:
522f6dfbefdSBarry Smith   This function should be passed to `SNESLineSearchSetPreCheck()`
52386d74e61SPeter Brune 
52486d74e61SPeter Brune   The justification for this method involves the linear convergence of a Picard iteration
5251d27aa22SBarry Smith   so the Picard linearization should be provided in place of the "Jacobian"  {cite}`hindmarsh1996time`. This correction
52686d74e61SPeter Brune   is generally not useful when using a Newton linearization.
52786d74e61SPeter Brune 
528420bcc1bSBarry Smith   Developer Note:
529420bcc1bSBarry Smith   The use of `PetscObjectGetState()` would eliminate the need for the `changed` argument to be provided
530420bcc1bSBarry Smith 
531420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESSetPicard()`, `SNESGetLineSearch()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchSetPostCheck()`
53286d74e61SPeter Brune @*/
533d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed, void *ctx)
534d71ae5a4SJacob Faibussowitsch {
53586d74e61SPeter Brune   PetscReal   angle = *(PetscReal *)linesearch->precheckctx;
53686d74e61SPeter Brune   Vec         Ylast;
53786d74e61SPeter Brune   PetscScalar dot;
53886d74e61SPeter Brune   PetscInt    iter;
53986d74e61SPeter Brune   PetscReal   ynorm, ylastnorm, theta, angle_radians;
54086d74e61SPeter Brune   SNES        snes;
54186d74e61SPeter Brune 
54286d74e61SPeter Brune   PetscFunctionBegin;
5439566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
5449566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)snes, "SNESLineSearchPreCheckPicard_Ylast", (PetscObject *)&Ylast));
54586d74e61SPeter Brune   if (!Ylast) {
5469566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(Y, &Ylast));
5479566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)snes, "SNESLineSearchPreCheckPicard_Ylast", (PetscObject)Ylast));
5489566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)Ylast));
54986d74e61SPeter Brune   }
5509566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &iter));
55186d74e61SPeter Brune   if (iter < 2) {
5529566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y, Ylast));
55386d74e61SPeter Brune     *changed = PETSC_FALSE;
5543ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
55586d74e61SPeter Brune   }
55686d74e61SPeter Brune 
5579566063dSJacob Faibussowitsch   PetscCall(VecDot(Y, Ylast, &dot));
5589566063dSJacob Faibussowitsch   PetscCall(VecNorm(Y, NORM_2, &ynorm));
5599566063dSJacob Faibussowitsch   PetscCall(VecNorm(Ylast, NORM_2, &ylastnorm));
560c197fec0SMatthew Knepley   if (ynorm == 0. || ylastnorm == 0.) {
561c197fec0SMatthew Knepley     *changed = PETSC_FALSE;
562c197fec0SMatthew Knepley     PetscFunctionReturn(PETSC_SUCCESS);
563c197fec0SMatthew Knepley   }
56486d74e61SPeter Brune   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
565255453a1SBarry Smith   theta         = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm), -1.0, 1.0));
56686d74e61SPeter Brune   angle_radians = angle * PETSC_PI / 180.;
56786d74e61SPeter Brune   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
56886d74e61SPeter Brune     /* Modify the step Y */
56986d74e61SPeter Brune     PetscReal alpha, ydiffnorm;
5709566063dSJacob Faibussowitsch     PetscCall(VecAXPY(Ylast, -1.0, Y));
5719566063dSJacob Faibussowitsch     PetscCall(VecNorm(Ylast, NORM_2, &ydiffnorm));
572f85e2ce2SBarry Smith     alpha = (ydiffnorm > .001 * ylastnorm) ? ylastnorm / ydiffnorm : 1000.0;
5739566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y, Ylast));
5749566063dSJacob Faibussowitsch     PetscCall(VecScale(Y, alpha));
575835f2295SStefano Zampini     PetscCall(PetscInfo(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));
576a47ec85fSBarry Smith     *changed = PETSC_TRUE;
57786d74e61SPeter Brune   } else {
578835f2295SStefano Zampini     PetscCall(PetscInfo(snes, "Angle %14.12e degrees exceeds threshold %14.12e, no correction applied\n", (double)(theta * 180 / PETSC_PI), (double)angle));
5799566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y, Ylast));
58086d74e61SPeter Brune     *changed = PETSC_FALSE;
581bf7f4e0aSPeter Brune   }
5823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
583bf7f4e0aSPeter Brune }
584bf7f4e0aSPeter Brune 
585f40b411bSPeter Brune /*@
586cd7522eaSPeter Brune   SNESLineSearchApply - Computes the line-search update.
587f40b411bSPeter Brune 
588c3339decSBarry Smith   Collective
589f40b411bSPeter Brune 
5906b095a85SStefano Zampini   Input Parameter:
5916b095a85SStefano Zampini . linesearch - The line search context
592f40b411bSPeter Brune 
5936b867d5aSJose E. Roman   Input/Output Parameters:
5946b867d5aSJose E. Roman + X     - The current solution, on output the new solution
595420bcc1bSBarry Smith . F     - The current function value, on output the new function value at the solution value `X`
5961bd63e3eSSatish Balay . fnorm - The current norm of `F`, on output the new norm of `F`
597a99ef635SJonas Heinzmann - Y     - The current search direction, on output the direction determined by the linesearch, i.e. `Xnew = Xold - lambda*Y`
598f40b411bSPeter Brune 
599cd7522eaSPeter Brune   Options Database Keys:
600a99ef635SJonas Heinzmann + -snes_linesearch_type                - basic (or equivalently none), bt, secant, cp, nleqerr, bisection, shell
601dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches
6021fe24845SBarry Smith . -snes_linesearch_damping             - The linesearch damping parameter, default is 1.0 (no damping)
6031fe24845SBarry Smith . -snes_linesearch_norms               - Turn on/off the linesearch norms computation (SNESLineSearchSetComputeNorms())
604a99ef635SJonas Heinzmann . -snes_linesearch_keeplambda          - Keep the previous `lambda` as the initial guess
6053c7d6663SPeter Brune - -snes_linesearch_max_it              - The number of iterations for iterative line searches
606cd7522eaSPeter Brune 
607e4094ef1SJacob Faibussowitsch   Level: intermediate
60820f4b53cSBarry Smith 
609cd7522eaSPeter Brune   Notes:
610f6dfbefdSBarry Smith   This is typically called from within a `SNESSolve()` implementation in order to
611f6dfbefdSBarry Smith   help with convergence of the nonlinear method.  Various `SNES` types use line searches
612cd7522eaSPeter Brune   in different ways, but the overarching theme is that a line search is used to determine
613a99ef635SJonas Heinzmann   an optimal damping parameter (that is `lambda`) of a step at each iteration of the method. Each
614f6dfbefdSBarry Smith   application of the line search may invoke `SNESComputeFunction()` several times, and
615cd7522eaSPeter Brune   therefore may be fairly expensive.
616cd7522eaSPeter Brune 
6171bd63e3eSSatish Balay .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchGetLambda()`, `SNESLineSearchPreCheck()`, `SNESLineSearchPostCheck()`, `SNESSolve()`, `SNESComputeFunction()`, `SNESLineSearchSetComputeNorms()`,
618db781477SPatrick Sanan           `SNESLineSearchType`, `SNESLineSearchSetType()`
619f40b411bSPeter Brune @*/
620d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal *fnorm, Vec Y)
621d71ae5a4SJacob Faibussowitsch {
622bf388a1fSBarry Smith   PetscFunctionBegin;
623f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
624bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
625bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
626064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(Y, VEC_CLASSID, 5);
627bf7f4e0aSPeter Brune 
628422a814eSBarry Smith   linesearch->result = SNES_LINESEARCH_SUCCEEDED;
629bf7f4e0aSPeter Brune 
630bf7f4e0aSPeter Brune   linesearch->vec_sol    = X;
631bf7f4e0aSPeter Brune   linesearch->vec_update = Y;
632bf7f4e0aSPeter Brune   linesearch->vec_func   = F;
633bf7f4e0aSPeter Brune 
6349566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchSetUp(linesearch));
635bf7f4e0aSPeter Brune 
636f5af7f23SKarl Rupp   if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */
637bf7f4e0aSPeter Brune 
638f5af7f23SKarl Rupp   if (fnorm) linesearch->fnorm = *fnorm;
63948a46eb9SPierre Jolivet   else PetscCall(VecNorm(F, NORM_2, &linesearch->fnorm));
640bf7f4e0aSPeter Brune 
6419566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(SNESLINESEARCH_Apply, linesearch, X, F, Y));
642bf7f4e0aSPeter Brune 
643dbbe0bcdSBarry Smith   PetscUseTypeMethod(linesearch, apply);
644bf7f4e0aSPeter Brune 
6459566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(SNESLINESEARCH_Apply, linesearch, X, F, Y));
646bf7f4e0aSPeter Brune 
647f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
6483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
649bf7f4e0aSPeter Brune }
650bf7f4e0aSPeter Brune 
651f40b411bSPeter Brune /*@
652f1c6b773SPeter Brune   SNESLineSearchDestroy - Destroys the line search instance.
653f40b411bSPeter Brune 
654c3339decSBarry Smith   Collective
655f40b411bSPeter Brune 
6562fe279fdSBarry Smith   Input Parameter:
6578e557f58SPeter Brune . linesearch - The line search context
658f40b411bSPeter Brune 
65984238204SBarry Smith   Level: developer
660f40b411bSPeter Brune 
661420bcc1bSBarry Smith   Note:
662420bcc1bSBarry Smith   The line search in `SNES` is automatically called on `SNESDestroy()` so this call is rarely needed
663420bcc1bSBarry Smith 
664420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchReset()`, `SNESDestroy()`
665f40b411bSPeter Brune @*/
666d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchDestroy(SNESLineSearch *linesearch)
667d71ae5a4SJacob Faibussowitsch {
668bf7f4e0aSPeter Brune   PetscFunctionBegin;
6693ba16761SJacob Faibussowitsch   if (!*linesearch) PetscFunctionReturn(PETSC_SUCCESS);
670f4f49eeaSPierre Jolivet   PetscValidHeaderSpecific(*linesearch, SNESLINESEARCH_CLASSID, 1);
671f4f49eeaSPierre Jolivet   if (--((PetscObject)*linesearch)->refct > 0) {
6729371c9d4SSatish Balay     *linesearch = NULL;
6733ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
6749371c9d4SSatish Balay   }
6759566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsViewOff((PetscObject)*linesearch));
6769566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchReset(*linesearch));
6773ba16761SJacob Faibussowitsch   PetscTryTypeMethod(*linesearch, destroy);
678648c30bcSBarry Smith   PetscCall(PetscViewerDestroy(&(*linesearch)->monitor));
679f4f49eeaSPierre Jolivet   PetscCall(SNESLineSearchMonitorCancel(*linesearch));
6809566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(linesearch));
6813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
682bf7f4e0aSPeter Brune }
683bf7f4e0aSPeter Brune 
684f40b411bSPeter Brune /*@
685dcf2fd19SBarry Smith   SNESLineSearchSetDefaultMonitor - Turns on/off printing useful information and debugging output about the line search.
686bf7f4e0aSPeter Brune 
687c3339decSBarry Smith   Logically Collective
688f6dfbefdSBarry Smith 
689bf7f4e0aSPeter Brune   Input Parameters:
690dcf2fd19SBarry Smith + linesearch - the linesearch object
69120f4b53cSBarry Smith - viewer     - an `PETSCVIEWERASCII` `PetscViewer` or `NULL` to turn off monitor
692bf7f4e0aSPeter Brune 
693f6dfbefdSBarry Smith   Options Database Key:
694dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - enables the monitor
695bf7f4e0aSPeter Brune 
696bf7f4e0aSPeter Brune   Level: intermediate
697bf7f4e0aSPeter Brune 
698e4094ef1SJacob Faibussowitsch   Developer Notes:
699f6dfbefdSBarry Smith   This monitor is implemented differently than the other line search monitors that are set with
700f6dfbefdSBarry Smith   `SNESLineSearchMonitorSet()` since it is called in many locations of the line search routines to display aspects of the
701d12e167eSBarry Smith   line search that are not visible to the other monitors.
702dcf2fd19SBarry Smith 
703420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `PETSCVIEWERASCII`, `SNESGetLineSearch()`, `SNESLineSearchGetDefaultMonitor()`, `PetscViewer`, `SNESLineSearchSetMonitor()`,
704f6dfbefdSBarry Smith           `SNESLineSearchMonitorSetFromOptions()`
705bf7f4e0aSPeter Brune @*/
706d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetDefaultMonitor(SNESLineSearch linesearch, PetscViewer viewer)
707d71ae5a4SJacob Faibussowitsch {
708bf7f4e0aSPeter Brune   PetscFunctionBegin;
709648c30bcSBarry Smith   PetscCall(PetscViewerDestroy(&linesearch->monitor));
710dcf2fd19SBarry Smith   linesearch->monitor = viewer;
7113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
712bf7f4e0aSPeter Brune }
713bf7f4e0aSPeter Brune 
714f40b411bSPeter Brune /*@
715420bcc1bSBarry Smith   SNESLineSearchGetDefaultMonitor - Gets the `PetscViewer` instance for the default line search monitor that is turned on with `SNESLineSearchSetDefaultMonitor()`
716f6dfbefdSBarry Smith 
717c3339decSBarry Smith   Logically Collective
7186a388c36SPeter Brune 
719f190f2fcSBarry Smith   Input Parameter:
720420bcc1bSBarry Smith . linesearch - the line search context
721f40b411bSPeter Brune 
722f190f2fcSBarry Smith   Output Parameter:
7238e557f58SPeter Brune . monitor - monitor context
724f40b411bSPeter Brune 
725f40b411bSPeter Brune   Level: intermediate
726f40b411bSPeter Brune 
727420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchSetDefaultMonitor()`, `PetscViewer`
728f40b411bSPeter Brune @*/
729d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetDefaultMonitor(SNESLineSearch linesearch, PetscViewer *monitor)
730d71ae5a4SJacob Faibussowitsch {
7316a388c36SPeter Brune   PetscFunctionBegin;
732f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
7336a388c36SPeter Brune   *monitor = linesearch->monitor;
7343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7356a388c36SPeter Brune }
7366a388c36SPeter Brune 
737dcf2fd19SBarry Smith /*@C
738420bcc1bSBarry Smith   SNESLineSearchMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated in the options database
739dcf2fd19SBarry Smith 
740c3339decSBarry Smith   Collective
741dcf2fd19SBarry Smith 
742dcf2fd19SBarry Smith   Input Parameters:
743420bcc1bSBarry Smith + ls           - `SNESLineSearch` object to monitor
744420bcc1bSBarry Smith . name         - the monitor type
745dcf2fd19SBarry Smith . help         - message indicating what monitoring is done
746dcf2fd19SBarry Smith . manual       - manual page for the monitor
74749abdd8aSBarry Smith . monitor      - the monitor function, must use `PetscViewerAndFormat` as its context
748f6dfbefdSBarry Smith - monitorsetup - a function that is called once ONLY if the user selected this monitor that may set additional features of the `SNESLineSearch` or `PetscViewer`
749dcf2fd19SBarry Smith 
750420bcc1bSBarry Smith   Calling sequence of `monitor`:
751420bcc1bSBarry Smith + ls - `SNESLineSearch` object being monitored
752420bcc1bSBarry Smith - vf - a `PetscViewerAndFormat` struct that provides the `PetscViewer` and `PetscViewerFormat` being used
753dcf2fd19SBarry Smith 
754420bcc1bSBarry Smith   Calling sequence of `monitorsetup`:
755420bcc1bSBarry Smith + ls - `SNESLineSearch` object being monitored
756420bcc1bSBarry Smith - vf - a `PetscViewerAndFormat` struct that provides the `PetscViewer` and `PetscViewerFormat` being used
757420bcc1bSBarry Smith 
758420bcc1bSBarry Smith   Level: advanced
759420bcc1bSBarry Smith 
760648c30bcSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetMonitor()`, `PetscOptionsCreateViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
761db781477SPatrick Sanan           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
762e4094ef1SJacob Faibussowitsch           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
763db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
764c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
765db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
766db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
767dcf2fd19SBarry Smith @*/
768420bcc1bSBarry Smith PetscErrorCode SNESLineSearchMonitorSetFromOptions(SNESLineSearch ls, const char name[], const char help[], const char manual[], PetscErrorCode (*monitor)(SNESLineSearch ls, PetscViewerAndFormat *vf), PetscErrorCode (*monitorsetup)(SNESLineSearch ls, PetscViewerAndFormat *vf))
769d71ae5a4SJacob Faibussowitsch {
770dcf2fd19SBarry Smith   PetscViewer       viewer;
771dcf2fd19SBarry Smith   PetscViewerFormat format;
772dcf2fd19SBarry Smith   PetscBool         flg;
773dcf2fd19SBarry Smith 
774dcf2fd19SBarry Smith   PetscFunctionBegin;
775648c30bcSBarry Smith   PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)ls), ((PetscObject)ls)->options, ((PetscObject)ls)->prefix, name, &viewer, &format, &flg));
776dcf2fd19SBarry Smith   if (flg) {
777d12e167eSBarry Smith     PetscViewerAndFormat *vf;
7789566063dSJacob Faibussowitsch     PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf));
779648c30bcSBarry Smith     PetscCall(PetscViewerDestroy(&viewer));
7801baa6e33SBarry Smith     if (monitorsetup) PetscCall((*monitorsetup)(ls, vf));
78149abdd8aSBarry Smith     PetscCall(SNESLineSearchMonitorSet(ls, (PetscErrorCode (*)(SNESLineSearch, void *))monitor, vf, (PetscCtxDestroyFn *)PetscViewerAndFormatDestroy));
782dcf2fd19SBarry Smith   }
7833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
784dcf2fd19SBarry Smith }
785dcf2fd19SBarry Smith 
786f40b411bSPeter Brune /*@
787f1c6b773SPeter Brune   SNESLineSearchSetFromOptions - Sets options for the line search
788f40b411bSPeter Brune 
789c3339decSBarry Smith   Logically Collective
790f6dfbefdSBarry Smith 
791f899ff85SJose E. Roman   Input Parameter:
792420bcc1bSBarry Smith . linesearch - a `SNESLineSearch` line search context
793f40b411bSPeter Brune 
794f40b411bSPeter Brune   Options Database Keys:
795a99ef635SJonas Heinzmann + -snes_linesearch_type <type>                                      - basic (or equivalently none), `bt`, `secant`, `cp`, `nleqerr`, `bisection`, `shell`
796a99ef635SJonas Heinzmann . -snes_linesearch_order <order>                                    - 1, 2, 3.  Most types only support certain orders (`bt` supports 1, 2 or 3)
797f6dfbefdSBarry Smith . -snes_linesearch_norms                                            - Turn on/off the linesearch norms for the basic linesearch typem (`SNESLineSearchSetComputeNorms()`)
798a99ef635SJonas Heinzmann . -snes_linesearch_minlambda                                        - The minimum `lambda`
799a99ef635SJonas Heinzmann . -snes_linesearch_maxlambda                                        - The maximum `lambda`
8001a9b3a06SPeter Brune . -snes_linesearch_rtol                                             - Relative tolerance for iterative line searches
8011a9b3a06SPeter Brune . -snes_linesearch_atol                                             - Absolute tolerance for iterative line searches
802a99ef635SJonas Heinzmann . -snes_linesearch_ltol                                             - Change in `lambda` tolerance for iterative line searches
8031a9b3a06SPeter Brune . -snes_linesearch_max_it                                           - The number of iterations for iterative line searches
804dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename]                              - Print progress of line searches
805dcf2fd19SBarry Smith . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine
8068e557f58SPeter Brune . -snes_linesearch_damping                                          - The linesearch damping parameter
807a99ef635SJonas Heinzmann . -snes_linesearch_keeplambda                                       - Keep the previous `lambda` as the initial guess.
8081a9b3a06SPeter Brune . -snes_linesearch_precheck_picard                                  - Use precheck that speeds up convergence of picard method
809d8d34be6SBarry Smith - -snes_linesearch_precheck_picard_angle                            - Angle used in Picard precheck method
810f40b411bSPeter Brune 
811f40b411bSPeter Brune   Level: intermediate
812f40b411bSPeter Brune 
81362842358SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchSetOrder()`, `SNESLineSearchSetType()`, `SNESLineSearchSetTolerances()`, `SNESLineSearchSetDamping()`, `SNESLineSearchPreCheckPicard()`,
814db781477SPatrick Sanan           `SNESLineSearchType`, `SNESLineSearchSetComputeNorms()`
815f40b411bSPeter Brune @*/
816d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch)
817d71ae5a4SJacob Faibussowitsch {
8181a4f838cSPeter Brune   const char *deft = SNESLINESEARCHBASIC;
819bf7f4e0aSPeter Brune   char        type[256];
820bf7f4e0aSPeter Brune   PetscBool   flg, set;
821dcf2fd19SBarry Smith   PetscViewer viewer;
822bf388a1fSBarry Smith 
823bf7f4e0aSPeter Brune   PetscFunctionBegin;
8249566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchRegisterAll());
825bf7f4e0aSPeter Brune 
826d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)linesearch);
827f5af7f23SKarl Rupp   if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name;
8289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-snes_linesearch_type", "Linesearch type", "SNESLineSearchSetType", SNESLineSearchList, deft, type, 256, &flg));
829bf7f4e0aSPeter Brune   if (flg) {
8309566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetType(linesearch, type));
831bf7f4e0aSPeter Brune   } else if (!((PetscObject)linesearch)->type_name) {
8329566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetType(linesearch, deft));
833bf7f4e0aSPeter Brune   }
834bf7f4e0aSPeter Brune 
835648c30bcSBarry Smith   PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)linesearch), ((PetscObject)linesearch)->options, ((PetscObject)linesearch)->prefix, "-snes_linesearch_monitor", &viewer, NULL, &set));
836648c30bcSBarry Smith   if (set) PetscCall(SNESLineSearchSetDefaultMonitor(linesearch, viewer));
8379566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchMonitorSetFromOptions(linesearch, "-snes_linesearch_monitor_solution_update", "View correction at each iteration", "SNESLineSearchMonitorSolutionUpdate", SNESLineSearchMonitorSolutionUpdate, NULL));
838bf7f4e0aSPeter Brune 
8391a9b3a06SPeter Brune   /* tolerances */
840a99ef635SJonas Heinzmann   PetscCall(PetscOptionsReal("-snes_linesearch_minlambda", "Minimum lambda", "SNESLineSearchSetTolerances", linesearch->minlambda, &linesearch->minlambda, NULL));
841a99ef635SJonas Heinzmann   PetscCall(PetscOptionsReal("-snes_linesearch_maxlambda", "Maximum lambda", "SNESLineSearchSetTolerances", linesearch->maxlambda, &linesearch->maxlambda, NULL));
8429566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_rtol", "Relative tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->rtol, &linesearch->rtol, NULL));
8439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_atol", "Absolute tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->atol, &linesearch->atol, NULL));
8449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_ltol", "Change in lambda tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->ltol, &linesearch->ltol, NULL));
845a99ef635SJonas Heinzmann   PetscCall(PetscOptionsInt("-snes_linesearch_max_it", "Maximum iterations for iterative line searches", "SNESLineSearchSetTolerances", linesearch->max_it, &linesearch->max_it, NULL));
846a99ef635SJonas Heinzmann 
847a99ef635SJonas Heinzmann   /* deprecated options */
848a99ef635SJonas Heinzmann   PetscCall(PetscOptionsDeprecated("-snes_linesearch_maxstep", "-snes_linesearch_maxlambda", "3.24.0", NULL));
8497a35526eSPeter Brune 
8501a9b3a06SPeter Brune   /* damping parameters */
851a99ef635SJonas Heinzmann   PetscCall(PetscOptionsReal("-snes_linesearch_damping", "Line search damping (and depending on chosen line search initial lambda guess)", "SNESLineSearchSetDamping", linesearch->damping, &linesearch->damping, NULL));
8521a9b3a06SPeter Brune 
8539566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_linesearch_keeplambda", "Use previous lambda as damping", "SNESLineSearchSetKeepLambda", linesearch->keeplambda, &linesearch->keeplambda, NULL));
8541a9b3a06SPeter Brune 
8551a9b3a06SPeter Brune   /* precheck */
8569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_linesearch_precheck_picard", "Use a correction that sometimes improves convergence of Picard iteration", "SNESLineSearchPreCheckPicard", flg, &flg, &set));
8577a35526eSPeter Brune   if (set) {
8587a35526eSPeter Brune     if (flg) {
8597a35526eSPeter Brune       linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */
860f5af7f23SKarl Rupp 
861d0609cedSBarry Smith       PetscCall(PetscOptionsReal("-snes_linesearch_precheck_picard_angle", "Maximum angle at which to activate the correction", "none", linesearch->precheck_picard_angle, &linesearch->precheck_picard_angle, NULL));
8629566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchSetPreCheck(linesearch, SNESLineSearchPreCheckPicard, &linesearch->precheck_picard_angle));
8637a35526eSPeter Brune     } else {
8649566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchSetPreCheck(linesearch, NULL, NULL));
8657a35526eSPeter Brune     }
8667a35526eSPeter Brune   }
8679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_linesearch_order", "Order of approximation used in the line search", "SNESLineSearchSetOrder", linesearch->order, &linesearch->order, NULL));
8689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_linesearch_norms", "Compute final norms in line search", "SNESLineSearchSetComputeNorms", linesearch->norms, &linesearch->norms, NULL));
8697a35526eSPeter Brune 
870dbbe0bcdSBarry Smith   PetscTryTypeMethod(linesearch, setfromoptions, PetscOptionsObject);
8715a9b6599SPeter Brune 
872dbbe0bcdSBarry Smith   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)linesearch, PetscOptionsObject));
873d0609cedSBarry Smith   PetscOptionsEnd();
8743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
875bf7f4e0aSPeter Brune }
876bf7f4e0aSPeter Brune 
877f40b411bSPeter Brune /*@
878f190f2fcSBarry Smith   SNESLineSearchView - Prints useful information about the line search
879f40b411bSPeter Brune 
88020f4b53cSBarry Smith   Logically Collective
88120f4b53cSBarry Smith 
882f40b411bSPeter Brune   Input Parameters:
8832fe279fdSBarry Smith + linesearch - line search context
884420bcc1bSBarry Smith - viewer     - the `PetscViewer` to display the line search information to
885f40b411bSPeter Brune 
886f40b411bSPeter Brune   Level: intermediate
887f40b411bSPeter Brune 
888420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `PetscViewer`, `SNESLineSearchCreate()`
889f40b411bSPeter Brune @*/
890d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer)
891d71ae5a4SJacob Faibussowitsch {
8929f196a02SMartin Diehl   PetscBool isascii;
893bf388a1fSBarry Smith 
894bf7f4e0aSPeter Brune   PetscFunctionBegin;
8957f1410a3SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
89648a46eb9SPierre Jolivet   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch), &viewer));
8977f1410a3SPeter Brune   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
8987f1410a3SPeter Brune   PetscCheckSameComm(linesearch, 1, viewer, 2);
899f40b411bSPeter Brune 
9009f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
9019f196a02SMartin Diehl   if (isascii) {
9029566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)linesearch, viewer));
9039566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
904dbbe0bcdSBarry Smith     PetscTryTypeMethod(linesearch, view, viewer);
9059566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
906a99ef635SJonas Heinzmann     PetscCall(PetscViewerASCIIPrintf(viewer, "  maxlambda=%e, minlambda=%e\n", (double)linesearch->maxlambda, (double)linesearch->minlambda));
9079566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol, (double)linesearch->atol, (double)linesearch->ltol));
908a99ef635SJonas Heinzmann     PetscCall(PetscViewerASCIIPrintf(viewer, "  maximum iterations=%" PetscInt_FMT "\n", linesearch->max_it));
9096b2b7091SBarry Smith     if (linesearch->ops->precheck) {
9106b2b7091SBarry Smith       if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) {
91163a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  using precheck step to speed up Picard convergence\n"));
9127f1410a3SPeter Brune       } else {
91363a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  using user-defined precheck step\n"));
9147f1410a3SPeter Brune       }
9157f1410a3SPeter Brune     }
91648a46eb9SPierre Jolivet     if (linesearch->ops->postcheck) PetscCall(PetscViewerASCIIPrintf(viewer, "  using user-defined postcheck step\n"));
9177f1410a3SPeter Brune   }
9183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
919bf7f4e0aSPeter Brune }
920bf7f4e0aSPeter Brune 
921cc4c1da9SBarry Smith /*@
922420bcc1bSBarry Smith   SNESLineSearchGetType - Gets the `SNESLinesearchType` of a `SNESLineSearch`
923a80ff896SJed Brown 
924c3339decSBarry Smith   Logically Collective
925a80ff896SJed Brown 
9262fe279fdSBarry Smith   Input Parameter:
927420bcc1bSBarry Smith . linesearch - the line search context
928a80ff896SJed Brown 
9292fe279fdSBarry Smith   Output Parameter:
9302fe279fdSBarry Smith . type - The type of line search, or `NULL` if not set
931a80ff896SJed Brown 
932a80ff896SJed Brown   Level: intermediate
933a80ff896SJed Brown 
934420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchSetType()`
935a80ff896SJed Brown @*/
936d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetType(SNESLineSearch linesearch, SNESLineSearchType *type)
937d71ae5a4SJacob Faibussowitsch {
938a80ff896SJed Brown   PetscFunctionBegin;
939a80ff896SJed Brown   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
9404f572ea9SToby Isaac   PetscAssertPointer(type, 2);
941a80ff896SJed Brown   *type = ((PetscObject)linesearch)->type_name;
9423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
943a80ff896SJed Brown }
944a80ff896SJed Brown 
945cc4c1da9SBarry Smith /*@
9460b4b7b1cSBarry Smith   SNESLineSearchSetType - Sets the `SNESLinesearchType` of a `SNESLineSearch` object to indicate the line search algorithm that should be used by a given `SNES` solver
947f40b411bSPeter Brune 
948c3339decSBarry Smith   Logically Collective
949f190f2fcSBarry Smith 
950f40b411bSPeter Brune   Input Parameters:
951420bcc1bSBarry Smith + linesearch - the line search context
952ceaaa498SBarry Smith - type       - The type of line search to be used, see `SNESLineSearchType`
9531fe24845SBarry Smith 
9543c7db156SBarry Smith   Options Database Key:
955a99ef635SJonas Heinzmann . -snes_linesearch_type <type> - basic (or equivalently none), bt, secant, cp, nleqerr, bisection, shell
956cd7522eaSPeter Brune 
957f40b411bSPeter Brune   Level: intermediate
958f40b411bSPeter Brune 
9590b4b7b1cSBarry Smith   Note:
9600b4b7b1cSBarry Smith   The `SNESLineSearch` object is generally obtained with `SNESGetLineSearch()`
9610b4b7b1cSBarry Smith 
9620b4b7b1cSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchGetType()`,
9630b4b7b1cSBarry Smith           `SNESGetLineSearch()`
964f40b411bSPeter Brune @*/
965d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type)
966d71ae5a4SJacob Faibussowitsch {
967bf7f4e0aSPeter Brune   PetscBool match;
9685f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(SNESLineSearch);
969bf7f4e0aSPeter Brune 
970bf7f4e0aSPeter Brune   PetscFunctionBegin;
971f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
9724f572ea9SToby Isaac   PetscAssertPointer(type, 2);
973bf7f4e0aSPeter Brune 
9749566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)linesearch, type, &match));
9753ba16761SJacob Faibussowitsch   if (match) PetscFunctionReturn(PETSC_SUCCESS);
976bf7f4e0aSPeter Brune 
9779566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(SNESLineSearchList, type, &r));
9786adde796SStefano Zampini   PetscCheck(r, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested Line Search type %s", type);
979bf7f4e0aSPeter Brune   /* Destroy the previous private line search context */
980213acdd3SPierre Jolivet   PetscTryTypeMethod(linesearch, destroy);
9810298fd71SBarry Smith   linesearch->ops->destroy = NULL;
982f1c6b773SPeter Brune   /* Reinitialize function pointers in SNESLineSearchOps structure */
9839e5d0892SLisandro Dalcin   linesearch->ops->apply          = NULL;
9849e5d0892SLisandro Dalcin   linesearch->ops->view           = NULL;
9859e5d0892SLisandro Dalcin   linesearch->ops->setfromoptions = NULL;
9869e5d0892SLisandro Dalcin   linesearch->ops->destroy        = NULL;
987bf7f4e0aSPeter Brune 
9889566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)linesearch, type));
9899566063dSJacob Faibussowitsch   PetscCall((*r)(linesearch));
9903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
991bf7f4e0aSPeter Brune }
992bf7f4e0aSPeter Brune 
993f40b411bSPeter Brune /*@
994f6dfbefdSBarry Smith   SNESLineSearchSetSNES - Sets the `SNES` for the linesearch for function evaluation.
995f40b411bSPeter Brune 
996f40b411bSPeter Brune   Input Parameters:
997420bcc1bSBarry Smith + linesearch - the line search context
998420bcc1bSBarry Smith - snes       - The `SNES` instance
999f40b411bSPeter Brune 
100078bcb3b5SPeter Brune   Level: developer
100178bcb3b5SPeter Brune 
1002f6dfbefdSBarry Smith   Note:
1003f190f2fcSBarry Smith   This happens automatically when the line search is obtained/created with
1004f6dfbefdSBarry Smith   `SNESGetLineSearch()`.  This routine is therefore mainly called within `SNES`
100578bcb3b5SPeter Brune   implementations.
1006f40b411bSPeter Brune 
1007420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetSNES()`, `SNESLineSearchSetVecs()`
1008f40b411bSPeter Brune @*/
1009d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes)
1010d71ae5a4SJacob Faibussowitsch {
1011bf7f4e0aSPeter Brune   PetscFunctionBegin;
1012f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1013bf7f4e0aSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 2);
1014bf7f4e0aSPeter Brune   linesearch->snes = snes;
10153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1016bf7f4e0aSPeter Brune }
1017bf7f4e0aSPeter Brune 
1018f40b411bSPeter Brune /*@
1019f6dfbefdSBarry Smith   SNESLineSearchGetSNES - Gets the `SNES` instance associated with the line search.
1020f6dfbefdSBarry Smith 
1021f6dfbefdSBarry Smith   Not Collective
1022f40b411bSPeter Brune 
10232fe279fdSBarry Smith   Input Parameter:
1024420bcc1bSBarry Smith . linesearch - the line search context
1025f40b411bSPeter Brune 
10262fe279fdSBarry Smith   Output Parameter:
10272fe279fdSBarry Smith . snes - The `SNES` instance
1028f40b411bSPeter Brune 
10298141a3b9SPeter Brune   Level: developer
1030f40b411bSPeter Brune 
1031420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESType`, `SNESLineSearchSetVecs()`
1032f40b411bSPeter Brune @*/
1033d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes)
1034d71ae5a4SJacob Faibussowitsch {
1035bf7f4e0aSPeter Brune   PetscFunctionBegin;
1036f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
10374f572ea9SToby Isaac   PetscAssertPointer(snes, 2);
1038bf7f4e0aSPeter Brune   *snes = linesearch->snes;
10393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1040bf7f4e0aSPeter Brune }
1041bf7f4e0aSPeter Brune 
1042f40b411bSPeter Brune /*@
1043a99ef635SJonas Heinzmann   SNESLineSearchGetLambda - Gets the last line search `lambda` used
1044f40b411bSPeter Brune 
1045f6dfbefdSBarry Smith   Not Collective
1046f6dfbefdSBarry Smith 
10472fe279fdSBarry Smith   Input Parameter:
1048420bcc1bSBarry Smith . linesearch - the line search context
1049f40b411bSPeter Brune 
10502fe279fdSBarry Smith   Output Parameter:
1051*8c5add6aSPierre Jolivet . lambda - The last `lambda` (scaling of the solution update) computed during `SNESLineSearchApply()`
1052f40b411bSPeter Brune 
105378bcb3b5SPeter Brune   Level: advanced
105478bcb3b5SPeter Brune 
1055f6dfbefdSBarry Smith   Note:
10568e557f58SPeter Brune   This is useful in methods where the solver is ill-scaled and
105778bcb3b5SPeter Brune   requires some adaptive notion of the difference in scale between the
1058f6dfbefdSBarry Smith   solution and the function.  For instance, `SNESQN` may be scaled by the
1059a99ef635SJonas Heinzmann   line search `lambda` using the argument -snes_qn_scaling ls.
106078bcb3b5SPeter Brune 
1061420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetLambda()`, `SNESLineSearchGetDamping()`, `SNESLineSearchApply()`
1062f40b411bSPeter Brune @*/
1063d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch linesearch, PetscReal *lambda)
1064d71ae5a4SJacob Faibussowitsch {
10656a388c36SPeter Brune   PetscFunctionBegin;
1066f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
10674f572ea9SToby Isaac   PetscAssertPointer(lambda, 2);
10686a388c36SPeter Brune   *lambda = linesearch->lambda;
10693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10706a388c36SPeter Brune }
10716a388c36SPeter Brune 
1072f40b411bSPeter Brune /*@
1073a99ef635SJonas Heinzmann   SNESLineSearchSetLambda - Sets the line search `lambda` (scaling of the solution update)
1074f40b411bSPeter Brune 
1075f40b411bSPeter Brune   Input Parameters:
10768e557f58SPeter Brune + linesearch - line search context
1077a99ef635SJonas Heinzmann - lambda     - The `lambda` to use
1078f40b411bSPeter Brune 
107920f4b53cSBarry Smith   Level: advanced
108020f4b53cSBarry Smith 
1081f6dfbefdSBarry Smith   Note:
1082f6dfbefdSBarry Smith   This routine is typically used within implementations of `SNESLineSearchApply()`
1083a99ef635SJonas Heinzmann   to set the final `lambda`.  This routine (and `SNESLineSearchGetLambda()`) were
1084a99ef635SJonas Heinzmann   added to facilitate Quasi-Newton methods that use the previous `lambda`
1085cd7522eaSPeter Brune   as an inner scaling parameter.
1086cd7522eaSPeter Brune 
1087420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetLambda()`
1088f40b411bSPeter Brune @*/
1089d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
1090d71ae5a4SJacob Faibussowitsch {
10916a388c36SPeter Brune   PetscFunctionBegin;
1092f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
10936a388c36SPeter Brune   linesearch->lambda = lambda;
10943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10956a388c36SPeter Brune }
10966a388c36SPeter Brune 
1097f40b411bSPeter Brune /*@
1098ceaaa498SBarry Smith   SNESLineSearchGetTolerances - Gets the tolerances for the line search.
1099f40b411bSPeter Brune 
1100f6dfbefdSBarry Smith   Not Collective
1101f6dfbefdSBarry Smith 
1102f899ff85SJose E. Roman   Input Parameter:
1103420bcc1bSBarry Smith . linesearch - the line search context
1104f40b411bSPeter Brune 
1105f40b411bSPeter Brune   Output Parameters:
1106a99ef635SJonas Heinzmann + minlambda - The minimum `lambda` allowed
1107a99ef635SJonas Heinzmann . maxlambda - The maximum `lambda` allowed
1108516fe3c3SPeter Brune . rtol      - The relative tolerance for iterative line searches
1109516fe3c3SPeter Brune . atol      - The absolute tolerance for iterative line searches
1110a99ef635SJonas Heinzmann . ltol      - The change in `lambda` tolerance for iterative line searches
1111a99ef635SJonas Heinzmann - max_it    - The maximum number of iterations of the line search
1112f40b411bSPeter Brune 
111378bcb3b5SPeter Brune   Level: intermediate
111478bcb3b5SPeter Brune 
1115f6dfbefdSBarry Smith   Note:
111678bcb3b5SPeter Brune   Different line searches may implement these parameters slightly differently as
11173c7d6663SPeter Brune   the type requires.
1118516fe3c3SPeter Brune 
1119420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetTolerances()`
1120f40b411bSPeter Brune @*/
1121a99ef635SJonas Heinzmann PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch linesearch, PetscReal *minlambda, PetscReal *maxlambda, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_it)
1122d71ae5a4SJacob Faibussowitsch {
11236a388c36SPeter Brune   PetscFunctionBegin;
1124f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1125a99ef635SJonas Heinzmann   if (minlambda) {
1126a99ef635SJonas Heinzmann     PetscAssertPointer(minlambda, 2);
1127a99ef635SJonas Heinzmann     *minlambda = linesearch->minlambda;
1128516fe3c3SPeter Brune   }
1129a99ef635SJonas Heinzmann   if (maxlambda) {
1130a99ef635SJonas Heinzmann     PetscAssertPointer(maxlambda, 3);
1131a99ef635SJonas Heinzmann     *maxlambda = linesearch->maxlambda;
1132516fe3c3SPeter Brune   }
1133516fe3c3SPeter Brune   if (rtol) {
11344f572ea9SToby Isaac     PetscAssertPointer(rtol, 4);
1135516fe3c3SPeter Brune     *rtol = linesearch->rtol;
1136516fe3c3SPeter Brune   }
1137516fe3c3SPeter Brune   if (atol) {
11384f572ea9SToby Isaac     PetscAssertPointer(atol, 5);
1139516fe3c3SPeter Brune     *atol = linesearch->atol;
1140516fe3c3SPeter Brune   }
1141516fe3c3SPeter Brune   if (ltol) {
11424f572ea9SToby Isaac     PetscAssertPointer(ltol, 6);
1143516fe3c3SPeter Brune     *ltol = linesearch->ltol;
1144516fe3c3SPeter Brune   }
1145a99ef635SJonas Heinzmann   if (max_it) {
1146a99ef635SJonas Heinzmann     PetscAssertPointer(max_it, 7);
1147a99ef635SJonas Heinzmann     *max_it = linesearch->max_it;
1148516fe3c3SPeter Brune   }
11493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11506a388c36SPeter Brune }
11516a388c36SPeter Brune 
1152f40b411bSPeter Brune /*@
1153ceaaa498SBarry Smith   SNESLineSearchSetTolerances -  Sets the tolerances for the linesearch.
1154f40b411bSPeter Brune 
1155f6dfbefdSBarry Smith   Collective
1156f6dfbefdSBarry Smith 
1157f40b411bSPeter Brune   Input Parameters:
1158420bcc1bSBarry Smith + linesearch - the line search context
1159a99ef635SJonas Heinzmann . minlambda  - The minimum `lambda` allowed
1160a99ef635SJonas Heinzmann . maxlambda  - The maximum `lambda` allowed
1161516fe3c3SPeter Brune . rtol       - The relative tolerance for iterative line searches
1162516fe3c3SPeter Brune . atol       - The absolute tolerance for iterative line searches
1163a99ef635SJonas Heinzmann . ltol       - The change in `lambda` tolerance for iterative line searches
1164420bcc1bSBarry Smith - max_it     - The maximum number of iterations of the line search
1165420bcc1bSBarry Smith 
1166420bcc1bSBarry Smith   Options Database Keys:
1167a99ef635SJonas Heinzmann + -snes_linesearch_minlambda - The minimum `lambda` allowed
1168a99ef635SJonas Heinzmann . -snes_linesearch_maxlambda - The maximum `lambda` allowed
1169420bcc1bSBarry Smith . -snes_linesearch_rtol      - Relative tolerance for iterative line searches
1170420bcc1bSBarry Smith . -snes_linesearch_atol      - Absolute tolerance for iterative line searches
1171a99ef635SJonas Heinzmann . -snes_linesearch_ltol      - Change in `lambda` tolerance for iterative line searches
1172420bcc1bSBarry Smith - -snes_linesearch_max_it    - The number of iterations for iterative line searches
1173f40b411bSPeter Brune 
117420f4b53cSBarry Smith   Level: intermediate
117520f4b53cSBarry Smith 
1176f6dfbefdSBarry Smith   Note:
1177420bcc1bSBarry Smith   The user may choose to not set any of the tolerances using `PETSC_DEFAULT` in place of an argument.
1178f40b411bSPeter Brune 
1179420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetTolerances()`
1180f40b411bSPeter Brune @*/
1181a99ef635SJonas Heinzmann PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch linesearch, PetscReal minlambda, PetscReal maxlambda, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_it)
1182d71ae5a4SJacob Faibussowitsch {
11836a388c36SPeter Brune   PetscFunctionBegin;
1184f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1185a99ef635SJonas Heinzmann   PetscValidLogicalCollectiveReal(linesearch, minlambda, 2);
1186a99ef635SJonas Heinzmann   PetscValidLogicalCollectiveReal(linesearch, maxlambda, 3);
1187d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, rtol, 4);
1188d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, atol, 5);
1189d3952184SSatish Balay   PetscValidLogicalCollectiveReal(linesearch, ltol, 6);
1190420bcc1bSBarry Smith   PetscValidLogicalCollectiveInt(linesearch, max_it, 7);
1191d3952184SSatish Balay 
1192a99ef635SJonas Heinzmann   if (minlambda != (PetscReal)PETSC_DEFAULT) {
1193a99ef635SJonas Heinzmann     PetscCheck(minlambda >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Minimum lambda %14.12e must be non-negative", (double)minlambda);
1194a99ef635SJonas Heinzmann     PetscCheck(minlambda < maxlambda, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Minimum lambda %14.12e must be smaller than maximum lambda %14.12e", (double)minlambda, (double)maxlambda);
1195a99ef635SJonas Heinzmann     linesearch->minlambda = minlambda;
1196d3952184SSatish Balay   }
1197d3952184SSatish Balay 
1198a99ef635SJonas Heinzmann   if (maxlambda != (PetscReal)PETSC_DEFAULT) {
1199a99ef635SJonas Heinzmann     PetscCheck(maxlambda > 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum lambda %14.12e must be positive", (double)maxlambda);
1200a99ef635SJonas Heinzmann     linesearch->maxlambda = maxlambda;
1201d3952184SSatish Balay   }
1202d3952184SSatish Balay 
120313bcc0bdSJacob Faibussowitsch   if (rtol != (PetscReal)PETSC_DEFAULT) {
12042061ca32SJunchao Zhang     PetscCheck(rtol >= 0.0 && rtol < 1.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Relative tolerance %14.12e must be non-negative and less than 1.0", (double)rtol);
1205516fe3c3SPeter Brune     linesearch->rtol = rtol;
1206d3952184SSatish Balay   }
1207d3952184SSatish Balay 
120813bcc0bdSJacob Faibussowitsch   if (atol != (PetscReal)PETSC_DEFAULT) {
12095f80ce2aSJacob Faibussowitsch     PetscCheck(atol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Absolute tolerance %14.12e must be non-negative", (double)atol);
1210516fe3c3SPeter Brune     linesearch->atol = atol;
1211d3952184SSatish Balay   }
1212d3952184SSatish Balay 
121313bcc0bdSJacob Faibussowitsch   if (ltol != (PetscReal)PETSC_DEFAULT) {
12145f80ce2aSJacob Faibussowitsch     PetscCheck(ltol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Lambda tolerance %14.12e must be non-negative", (double)ltol);
1215516fe3c3SPeter Brune     linesearch->ltol = ltol;
1216d3952184SSatish Balay   }
1217d3952184SSatish Balay 
1218420bcc1bSBarry Smith   if (max_it != PETSC_DEFAULT) {
1219420bcc1bSBarry Smith     PetscCheck(max_it >= 0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of iterations %" PetscInt_FMT " must be non-negative", max_it);
1220a99ef635SJonas Heinzmann     linesearch->max_it = max_it;
1221d3952184SSatish Balay   }
12223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12236a388c36SPeter Brune }
12246a388c36SPeter Brune 
1225f40b411bSPeter Brune /*@
1226f1c6b773SPeter Brune   SNESLineSearchGetDamping - Gets the line search damping parameter.
1227f40b411bSPeter Brune 
12282fe279fdSBarry Smith   Input Parameter:
1229420bcc1bSBarry Smith . linesearch - the line search context
1230f40b411bSPeter Brune 
12312fe279fdSBarry Smith   Output Parameter:
12328e557f58SPeter Brune . damping - The damping parameter
1233f40b411bSPeter Brune 
123478bcb3b5SPeter Brune   Level: advanced
1235f40b411bSPeter Brune 
1236420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearchGetStepTolerance()`, `SNESQN`
1237f40b411bSPeter Brune @*/
1238d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch linesearch, PetscReal *damping)
1239d71ae5a4SJacob Faibussowitsch {
12406a388c36SPeter Brune   PetscFunctionBegin;
1241f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
12424f572ea9SToby Isaac   PetscAssertPointer(damping, 2);
12436a388c36SPeter Brune   *damping = linesearch->damping;
12443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12456a388c36SPeter Brune }
12466a388c36SPeter Brune 
1247f40b411bSPeter Brune /*@
1248fd292e60Sprj-   SNESLineSearchSetDamping - Sets the line search damping parameter.
1249f40b411bSPeter Brune 
1250f40b411bSPeter Brune   Input Parameters:
1251420bcc1bSBarry Smith + linesearch - the line search context
125203fd4120SBarry Smith - damping    - The damping parameter
1253f40b411bSPeter Brune 
1254f6dfbefdSBarry Smith   Options Database Key:
1255f6dfbefdSBarry Smith . -snes_linesearch_damping <damping> - the damping value
1256f6dfbefdSBarry Smith 
1257f40b411bSPeter Brune   Level: intermediate
1258f40b411bSPeter Brune 
1259f6dfbefdSBarry Smith   Note:
1260f6dfbefdSBarry Smith   The `SNESLINESEARCHNONE` line search merely takes the update step scaled by the damping parameter.
1261a99ef635SJonas Heinzmann   The use of the damping parameter in the `SNESLINESEARCHSECANT` and `SNESLINESEARCHCP` line searches is much more subtle;
1262a99ef635SJonas Heinzmann   it is used as a starting point for the secant method. Depending on the choice for `maxlambda`,
1263a99ef635SJonas Heinzmann   the eventual `lambda` may be greater than the damping parameter however.
1264a99ef635SJonas Heinzmann   For `SNESLINESEARCHBISECTION` and `SNESLINESEARCHBT` the damping is instead used as the initial guess,
1265a99ef635SJonas Heinzmann   below which the line search will not go. Hence, it is the maximum possible value for `lambda`.
1266cd7522eaSPeter Brune 
1267420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetDamping()`
1268f40b411bSPeter Brune @*/
1269d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch linesearch, PetscReal damping)
1270d71ae5a4SJacob Faibussowitsch {
12716a388c36SPeter Brune   PetscFunctionBegin;
1272f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
12736a388c36SPeter Brune   linesearch->damping = damping;
12743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12756a388c36SPeter Brune }
12766a388c36SPeter Brune 
127759405d5eSPeter Brune /*@
127859405d5eSPeter Brune   SNESLineSearchGetOrder - Gets the line search approximation order.
127959405d5eSPeter Brune 
1280f6dfbefdSBarry Smith   Input Parameter:
1281420bcc1bSBarry Smith . linesearch - the line search context
128259405d5eSPeter Brune 
1283f6dfbefdSBarry Smith   Output Parameter:
12848e557f58SPeter Brune . order - The order
128559405d5eSPeter Brune 
128659405d5eSPeter Brune   Level: intermediate
128759405d5eSPeter Brune 
1288420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetOrder()`
128959405d5eSPeter Brune @*/
1290d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch linesearch, PetscInt *order)
1291d71ae5a4SJacob Faibussowitsch {
129259405d5eSPeter Brune   PetscFunctionBegin;
129359405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
12944f572ea9SToby Isaac   PetscAssertPointer(order, 2);
129559405d5eSPeter Brune   *order = linesearch->order;
12963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
129759405d5eSPeter Brune }
129859405d5eSPeter Brune 
129959405d5eSPeter Brune /*@
13001f8196a2SJed Brown   SNESLineSearchSetOrder - Sets the maximum order of the polynomial fit used in the line search
130159405d5eSPeter Brune 
130259405d5eSPeter Brune   Input Parameters:
1303420bcc1bSBarry Smith + linesearch - the line search context
1304ceaaa498SBarry Smith - order      - The order
130559405d5eSPeter Brune 
130659405d5eSPeter Brune   Level: intermediate
130759405d5eSPeter Brune 
1308ceaaa498SBarry Smith   Values for `order`\:
1309f6dfbefdSBarry Smith +  1 or `SNES_LINESEARCH_ORDER_LINEAR` - linear order
1310f6dfbefdSBarry Smith .  2 or `SNES_LINESEARCH_ORDER_QUADRATIC` - quadratic order
1311f6dfbefdSBarry Smith -  3 or `SNES_LINESEARCH_ORDER_CUBIC` - cubic order
131278bcb3b5SPeter Brune 
1313420bcc1bSBarry Smith   Options Database Key:
1314420bcc1bSBarry Smith . -snes_linesearch_order <order> - 1, 2, 3.  Most types only support certain orders (`SNESLINESEARCHBT` supports 2 or 3)
1315420bcc1bSBarry Smith 
1316ceaaa498SBarry Smith   Note:
1317ceaaa498SBarry Smith   These orders are supported by `SNESLINESEARCHBT` and `SNESLINESEARCHCP`
131859405d5eSPeter Brune 
1319420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetOrder()`, `SNESLineSearchSetDamping()`
132059405d5eSPeter Brune @*/
1321d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch linesearch, PetscInt order)
1322d71ae5a4SJacob Faibussowitsch {
132359405d5eSPeter Brune   PetscFunctionBegin;
132459405d5eSPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
132559405d5eSPeter Brune   linesearch->order = order;
13263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
132759405d5eSPeter Brune }
132859405d5eSPeter Brune 
1329f40b411bSPeter Brune /*@
1330420bcc1bSBarry Smith   SNESLineSearchGetNorms - Gets the norms for the current solution `X`, the current update `Y`, and the current function value `F`.
1331f40b411bSPeter Brune 
1332f6dfbefdSBarry Smith   Not Collective
1333f6dfbefdSBarry Smith 
1334f899ff85SJose E. Roman   Input Parameter:
1335420bcc1bSBarry Smith . linesearch - the line search context
1336f40b411bSPeter Brune 
1337f40b411bSPeter Brune   Output Parameters:
1338f40b411bSPeter Brune + xnorm - The norm of the current solution
1339a68bbae5SBarry Smith . fnorm - The norm of the current function, this is the `norm(function(X))` where `X` is the current solution.
1340a99ef635SJonas Heinzmann - ynorm - The norm of the current update (after scaling by the linesearch computed `lambda`)
1341f40b411bSPeter Brune 
134278bcb3b5SPeter Brune   Level: developer
1343f40b411bSPeter Brune 
1344420bcc1bSBarry Smith   Notes:
1345420bcc1bSBarry Smith   Some values may not be up-to-date at particular points in the code.
1346a68bbae5SBarry Smith 
1347a68bbae5SBarry Smith   This, in combination with `SNESLineSearchSetNorms()`, allow the line search and the `SNESSolve_XXX()` to share
1348a68bbae5SBarry Smith   computed values.
1349a68bbae5SBarry Smith 
13500241b274SPierre Jolivet .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetNorms()`, `SNESLineSearchGetVecs()`
1351f40b411bSPeter Brune @*/
1352d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal *xnorm, PetscReal *fnorm, PetscReal *ynorm)
1353d71ae5a4SJacob Faibussowitsch {
1354bf7f4e0aSPeter Brune   PetscFunctionBegin;
1355f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1356f5af7f23SKarl Rupp   if (xnorm) *xnorm = linesearch->xnorm;
1357f5af7f23SKarl Rupp   if (fnorm) *fnorm = linesearch->fnorm;
1358f5af7f23SKarl Rupp   if (ynorm) *ynorm = linesearch->ynorm;
13593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1360bf7f4e0aSPeter Brune }
1361bf7f4e0aSPeter Brune 
1362f40b411bSPeter Brune /*@
13631bd63e3eSSatish Balay   SNESLineSearchSetNorms - Sets the computed norms for the current solution `X`, the current update `Y`, and the current function value `F`.
1364f40b411bSPeter Brune 
1365c3339decSBarry Smith   Collective
1366f6dfbefdSBarry Smith 
1367f40b411bSPeter Brune   Input Parameters:
1368420bcc1bSBarry Smith + linesearch - the line search context
1369f40b411bSPeter Brune . xnorm      - The norm of the current solution
1370a68bbae5SBarry Smith . fnorm      - The norm of the current function, this is the `norm(function(X))` where `X` is the current solution
1371a99ef635SJonas Heinzmann - ynorm      - The norm of the current update (after scaling by the linesearch computed `lambda`)
1372f40b411bSPeter Brune 
1373f6dfbefdSBarry Smith   Level: developer
1374f40b411bSPeter Brune 
1375420bcc1bSBarry Smith   Note:
1376420bcc1bSBarry Smith   This is called by the line search routines to store the values they have just computed
1377420bcc1bSBarry Smith 
1378420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()`
1379f40b411bSPeter Brune @*/
1380d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
1381d71ae5a4SJacob Faibussowitsch {
13826a388c36SPeter Brune   PetscFunctionBegin;
1383f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
13846a388c36SPeter Brune   linesearch->xnorm = xnorm;
13856a388c36SPeter Brune   linesearch->fnorm = fnorm;
13866a388c36SPeter Brune   linesearch->ynorm = ynorm;
13873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13886a388c36SPeter Brune }
13896a388c36SPeter Brune 
1390f40b411bSPeter Brune /*@
1391420bcc1bSBarry Smith   SNESLineSearchComputeNorms - Explicitly computes the norms of the current solution `X`, the current update `Y`, and the current function value `F`.
1392f40b411bSPeter Brune 
13932fe279fdSBarry Smith   Input Parameter:
1394420bcc1bSBarry Smith . linesearch - the line search context
1395f40b411bSPeter Brune 
139620f4b53cSBarry Smith   Options Database Key:
13978e557f58SPeter Brune . -snes_linesearch_norms - turn norm computation on or off
1398f40b411bSPeter Brune 
1399f40b411bSPeter Brune   Level: intermediate
1400f40b411bSPeter Brune 
1401420bcc1bSBarry Smith   Developer Note:
1402420bcc1bSBarry Smith   The options database key is misnamed. It should be -snes_linesearch_compute_norms
1403420bcc1bSBarry Smith 
1404420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms`, `SNESLineSearchSetNorms()`, `SNESLineSearchSetComputeNorms()`
1405f40b411bSPeter Brune @*/
1406d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
1407d71ae5a4SJacob Faibussowitsch {
14089bd66eb0SPeter Brune   SNES snes;
1409bf388a1fSBarry Smith 
14106a388c36SPeter Brune   PetscFunctionBegin;
14116a388c36SPeter Brune   if (linesearch->norms) {
14129bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
14139566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
14149566063dSJacob Faibussowitsch       PetscCall(VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
14159566063dSJacob Faibussowitsch       PetscCall(VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm));
14169566063dSJacob Faibussowitsch       PetscCall((*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm));
14179bd66eb0SPeter Brune     } else {
14189566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(linesearch->vec_func, NORM_2, &linesearch->fnorm));
14199566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
14209566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm));
14219566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(linesearch->vec_func, NORM_2, &linesearch->fnorm));
14229566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
14239566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(linesearch->vec_update, NORM_2, &linesearch->ynorm));
14246a388c36SPeter Brune     }
14259bd66eb0SPeter Brune   }
14263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14276a388c36SPeter Brune }
14286a388c36SPeter Brune 
14296f263ca3SPeter Brune /*@
14306f263ca3SPeter Brune   SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search.
14316f263ca3SPeter Brune 
14326f263ca3SPeter Brune   Input Parameters:
1433420bcc1bSBarry Smith + linesearch - the line search context
143478bcb3b5SPeter Brune - flg        - indicates whether or not to compute norms
14356f263ca3SPeter Brune 
143620f4b53cSBarry Smith   Options Database Key:
1437420bcc1bSBarry Smith . -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic (none) `SNESLINESEARCHBASIC` line search
14386f263ca3SPeter Brune 
143920f4b53cSBarry Smith   Level: intermediate
144020f4b53cSBarry Smith 
1441f6dfbefdSBarry Smith   Note:
1442f6dfbefdSBarry Smith   This is most relevant to the `SNESLINESEARCHBASIC` (or equivalently `SNESLINESEARCHNONE`) line search type since most line searches have a stopping criteria involving the norm.
14436f263ca3SPeter Brune 
1444420bcc1bSBarry Smith   Developer Note:
1445420bcc1bSBarry Smith   The options database key is misnamed. It should be -snes_linesearch_compute_norms
1446420bcc1bSBarry Smith 
1447420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetNorms()`, `SNESLineSearchComputeNorms()`, `SNESLINESEARCHBASIC`
14486f263ca3SPeter Brune @*/
1449d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg)
1450d71ae5a4SJacob Faibussowitsch {
14516f263ca3SPeter Brune   PetscFunctionBegin;
14526f263ca3SPeter Brune   linesearch->norms = flg;
14533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14546f263ca3SPeter Brune }
14556f263ca3SPeter Brune 
1456f40b411bSPeter Brune /*@
1457f6dfbefdSBarry Smith   SNESLineSearchGetVecs - Gets the vectors from the `SNESLineSearch` context
1458f6dfbefdSBarry Smith 
14598f14a041SBarry Smith   Not Collective but the vectors are parallel
1460f40b411bSPeter Brune 
1461f899ff85SJose E. Roman   Input Parameter:
1462420bcc1bSBarry Smith . linesearch - the line search context
1463f40b411bSPeter Brune 
1464f40b411bSPeter Brune   Output Parameters:
14656232e825SPeter Brune + X - Solution vector
14666232e825SPeter Brune . F - Function vector
14676232e825SPeter Brune . Y - Search direction vector
14686232e825SPeter Brune . W - Solution work vector
14696232e825SPeter Brune - G - Function work vector
14706232e825SPeter Brune 
147120f4b53cSBarry Smith   Level: advanced
147220f4b53cSBarry Smith 
14737bba9028SPeter Brune   Notes:
147420f4b53cSBarry Smith   At the beginning of a line search application, `X` should contain a
147520f4b53cSBarry Smith   solution and the vector `F` the function computed at `X`.  At the end of the
147620f4b53cSBarry Smith   line search application, `X` should contain the new solution, and `F` the
14776232e825SPeter Brune   function evaluated at the new solution.
1478f40b411bSPeter Brune 
1479f6dfbefdSBarry Smith   These vectors are owned by the `SNESLineSearch` and should not be destroyed by the caller
14802a7a6963SBarry Smith 
1481420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()`
1482f40b411bSPeter Brune @*/
1483d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch, Vec *X, Vec *F, Vec *Y, Vec *W, Vec *G)
1484d71ae5a4SJacob Faibussowitsch {
14856a388c36SPeter Brune   PetscFunctionBegin;
1486f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
14876a388c36SPeter Brune   if (X) {
14884f572ea9SToby Isaac     PetscAssertPointer(X, 2);
14896a388c36SPeter Brune     *X = linesearch->vec_sol;
14906a388c36SPeter Brune   }
14916a388c36SPeter Brune   if (F) {
14924f572ea9SToby Isaac     PetscAssertPointer(F, 3);
14936a388c36SPeter Brune     *F = linesearch->vec_func;
14946a388c36SPeter Brune   }
14956a388c36SPeter Brune   if (Y) {
14964f572ea9SToby Isaac     PetscAssertPointer(Y, 4);
14976a388c36SPeter Brune     *Y = linesearch->vec_update;
14986a388c36SPeter Brune   }
14996a388c36SPeter Brune   if (W) {
15004f572ea9SToby Isaac     PetscAssertPointer(W, 5);
15016a388c36SPeter Brune     *W = linesearch->vec_sol_new;
15026a388c36SPeter Brune   }
15036a388c36SPeter Brune   if (G) {
15044f572ea9SToby Isaac     PetscAssertPointer(G, 6);
15056a388c36SPeter Brune     *G = linesearch->vec_func_new;
15066a388c36SPeter Brune   }
15073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15086a388c36SPeter Brune }
15096a388c36SPeter Brune 
1510f40b411bSPeter Brune /*@
1511f6dfbefdSBarry Smith   SNESLineSearchSetVecs - Sets the vectors on the `SNESLineSearch` context
1512f6dfbefdSBarry Smith 
1513c3339decSBarry Smith   Logically Collective
1514f40b411bSPeter Brune 
1515f40b411bSPeter Brune   Input Parameters:
1516420bcc1bSBarry Smith + linesearch - the line search context
15176232e825SPeter Brune . X          - Solution vector
15186232e825SPeter Brune . F          - Function vector
15196232e825SPeter Brune . Y          - Search direction vector
15206232e825SPeter Brune . W          - Solution work vector
15216232e825SPeter Brune - G          - Function work vector
1522f40b411bSPeter Brune 
1523420bcc1bSBarry Smith   Level: developer
1524f40b411bSPeter Brune 
1525420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetNorms()`, `SNESLineSearchGetVecs()`
1526f40b411bSPeter Brune @*/
1527d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch, Vec X, Vec F, Vec Y, Vec W, Vec G)
1528d71ae5a4SJacob Faibussowitsch {
15296a388c36SPeter Brune   PetscFunctionBegin;
1530f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
15316a388c36SPeter Brune   if (X) {
15326a388c36SPeter Brune     PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
15336a388c36SPeter Brune     linesearch->vec_sol = X;
15346a388c36SPeter Brune   }
15356a388c36SPeter Brune   if (F) {
15366a388c36SPeter Brune     PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
15376a388c36SPeter Brune     linesearch->vec_func = F;
15386a388c36SPeter Brune   }
15396a388c36SPeter Brune   if (Y) {
15406a388c36SPeter Brune     PetscValidHeaderSpecific(Y, VEC_CLASSID, 4);
15416a388c36SPeter Brune     linesearch->vec_update = Y;
15426a388c36SPeter Brune   }
15436a388c36SPeter Brune   if (W) {
15446a388c36SPeter Brune     PetscValidHeaderSpecific(W, VEC_CLASSID, 5);
15456a388c36SPeter Brune     linesearch->vec_sol_new = W;
15466a388c36SPeter Brune   }
15476a388c36SPeter Brune   if (G) {
15486a388c36SPeter Brune     PetscValidHeaderSpecific(G, VEC_CLASSID, 6);
15496a388c36SPeter Brune     linesearch->vec_func_new = G;
15506a388c36SPeter Brune   }
15513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15526a388c36SPeter Brune }
15536a388c36SPeter Brune 
1554cc4c1da9SBarry Smith /*@
1555f1c6b773SPeter Brune   SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1556f6dfbefdSBarry Smith   `SNESLineSearch` options in the database.
1557e7058c64SPeter Brune 
1558c3339decSBarry Smith   Logically Collective
1559e7058c64SPeter Brune 
1560e7058c64SPeter Brune   Input Parameters:
1561f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context
1562e7058c64SPeter Brune - prefix     - the prefix to prepend to all option names
1563e7058c64SPeter Brune 
156420f4b53cSBarry Smith   Level: advanced
156520f4b53cSBarry Smith 
1566f6dfbefdSBarry Smith   Note:
1567e7058c64SPeter Brune   A hyphen (-) must NOT be given at the beginning of the prefix name.
1568e7058c64SPeter Brune   The first character of all runtime options is AUTOMATICALLY the hyphen.
1569e7058c64SPeter Brune 
1570420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch()`, `SNESLineSearchSetFromOptions()`, `SNESGetOptionsPrefix()`
1571e7058c64SPeter Brune @*/
1572d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch, const char prefix[])
1573d71ae5a4SJacob Faibussowitsch {
1574e7058c64SPeter Brune   PetscFunctionBegin;
1575f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
15769566063dSJacob Faibussowitsch   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)linesearch, prefix));
15773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1578e7058c64SPeter Brune }
1579e7058c64SPeter Brune 
1580cc4c1da9SBarry Smith /*@
1581f6dfbefdSBarry Smith   SNESLineSearchGetOptionsPrefix - Gets the prefix used for searching for all
1582f1c6b773SPeter Brune   SNESLineSearch options in the database.
1583e7058c64SPeter Brune 
1584e7058c64SPeter Brune   Not Collective
1585e7058c64SPeter Brune 
1586e7058c64SPeter Brune   Input Parameter:
1587f6dfbefdSBarry Smith . linesearch - the `SNESLineSearch` context
1588e7058c64SPeter Brune 
1589e7058c64SPeter Brune   Output Parameter:
1590e7058c64SPeter Brune . prefix - pointer to the prefix string used
1591e7058c64SPeter Brune 
1592e7058c64SPeter Brune   Level: advanced
1593e7058c64SPeter Brune 
1594420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESAppendOptionsPrefix()`
1595e7058c64SPeter Brune @*/
1596d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch, const char *prefix[])
1597d71ae5a4SJacob Faibussowitsch {
1598e7058c64SPeter Brune   PetscFunctionBegin;
1599f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
16009566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)linesearch, prefix));
16013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1602e7058c64SPeter Brune }
1603bf7f4e0aSPeter Brune 
16048d359177SBarry Smith /*@C
1605f6dfbefdSBarry Smith   SNESLineSearchSetWorkVecs - Sets work vectors for the line search.
1606f40b411bSPeter Brune 
1607d8d19677SJose E. Roman   Input Parameters:
1608f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context
1609f40b411bSPeter Brune - nwork      - the number of work vectors
1610f40b411bSPeter Brune 
1611f40b411bSPeter Brune   Level: developer
1612f40b411bSPeter Brune 
1613420bcc1bSBarry Smith   Developer Note:
1614420bcc1bSBarry Smith   This is called from within the set up routines for each of the line search types `SNESLineSearchType`
1615420bcc1bSBarry Smith 
1616420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESSetWorkVecs()`
1617f40b411bSPeter Brune @*/
1618d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork)
1619d71ae5a4SJacob Faibussowitsch {
1620bf7f4e0aSPeter Brune   PetscFunctionBegin;
16210fdf79fbSJacob Faibussowitsch   PetscCheck(linesearch->vec_sol, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
16229566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work));
16233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1624bf7f4e0aSPeter Brune }
1625bf7f4e0aSPeter Brune 
1626f40b411bSPeter Brune /*@
1627422a814eSBarry Smith   SNESLineSearchGetReason - Gets the success/failure status of the last line search application
1628f40b411bSPeter Brune 
16292fe279fdSBarry Smith   Input Parameter:
1630420bcc1bSBarry Smith . linesearch - the line search context
1631f40b411bSPeter Brune 
16322fe279fdSBarry Smith   Output Parameter:
1633422a814eSBarry Smith . result - The success or failure status
1634f40b411bSPeter Brune 
163520f4b53cSBarry Smith   Level: developer
163620f4b53cSBarry Smith 
1637f6dfbefdSBarry Smith   Note:
1638420bcc1bSBarry Smith   This is typically called after `SNESLineSearchApply()` in order to determine if the line search failed
1639420bcc1bSBarry Smith   (and set into the `SNES` convergence accordingly).
1640cd7522eaSPeter Brune 
1641420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetReason()`, `SNESLineSearchReason`
1642f40b411bSPeter Brune @*/
1643d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result)
1644d71ae5a4SJacob Faibussowitsch {
1645bf7f4e0aSPeter Brune   PetscFunctionBegin;
1646f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
16474f572ea9SToby Isaac   PetscAssertPointer(result, 2);
1648422a814eSBarry Smith   *result = linesearch->result;
16493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1650bf7f4e0aSPeter Brune }
1651bf7f4e0aSPeter Brune 
1652f40b411bSPeter Brune /*@
1653420bcc1bSBarry Smith   SNESLineSearchSetReason - Sets the success/failure status of the line search application
1654f40b411bSPeter Brune 
16555d83a8b1SBarry Smith   Logically Collective; No Fortran Support
16565d83a8b1SBarry Smith 
1657f40b411bSPeter Brune   Input Parameters:
1658420bcc1bSBarry Smith + linesearch - the line search context
1659422a814eSBarry Smith - result     - The success or failure status
1660f40b411bSPeter Brune 
166120f4b53cSBarry Smith   Level: developer
166220f4b53cSBarry Smith 
1663f6dfbefdSBarry Smith   Note:
1664420bcc1bSBarry Smith   This is typically called in a `SNESLineSearchType` implementation of `SNESLineSearchApply()` or a `SNESLINESEARCHSHELL` implementation to set
1665cd7522eaSPeter Brune   the success or failure of the line search method.
1666cd7522eaSPeter Brune 
1667420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchReason`, `SNESLineSearchGetSResult()`
1668f40b411bSPeter Brune @*/
1669d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result)
1670d71ae5a4SJacob Faibussowitsch {
16716a388c36SPeter Brune   PetscFunctionBegin;
16725fd66863SKarl Rupp   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1673422a814eSBarry Smith   linesearch->result = result;
16743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16756a388c36SPeter Brune }
16766a388c36SPeter Brune 
1677ceaaa498SBarry Smith // PetscClangLinter pragma disable: -fdoc-param-list-func-parameter-documentation
16789bd66eb0SPeter Brune /*@C
1679f1c6b773SPeter Brune   SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.
16809bd66eb0SPeter Brune 
1681c3339decSBarry Smith   Logically Collective
1682f6dfbefdSBarry Smith 
16839bd66eb0SPeter Brune   Input Parameters:
1684ceaaa498SBarry Smith + linesearch   - the linesearch object
16858434afd1SBarry Smith . projectfunc  - function for projecting the function to the bounds, see `SNESLineSearchVIProjectFn` for calling sequence
1686d5def619SJonas Heinzmann . normfunc     - function for computing the norm of an active set, see `SNESLineSearchVINormFn` for calling sequence
1687d5def619SJonas Heinzmann - dirderivfunc - function for computing the directional derivative of an active set, see `SNESLineSearchVIDirDerivFn` for calling sequence
16889bd66eb0SPeter Brune 
1689f6dfbefdSBarry Smith   Level: advanced
16909bd66eb0SPeter Brune 
1691ceaaa498SBarry Smith   Notes:
169220f4b53cSBarry Smith   The VI solvers require projection of the solution to the feasible set.  `projectfunc` should implement this.
169320f4b53cSBarry Smith 
169420f4b53cSBarry Smith   The VI solvers require special evaluation of the function norm such that the norm is only calculated
169520f4b53cSBarry Smith   on the inactive set.  This should be implemented by `normfunc`.
169620f4b53cSBarry Smith 
1697d5def619SJonas Heinzmann   The VI solvers further require special evaluation of the directional derivative (when assuming that there exists some $G(x)$
1698d5def619SJonas Heinzmann   for which the `SNESFunctionFn` $F(x) = grad G(x)$) such that it is only calculated on the inactive set.
1699d5def619SJonas Heinzmann   This should be implemented by `dirderivfunc`.
1700d5def619SJonas Heinzmann 
17019bcc50f1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetVIFunctions()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()`,
1702d5def619SJonas Heinzmann           `SNESLineSearchVIProjectFn`, `SNESLineSearchVINormFn`, `SNESLineSearchVIDirDerivFn`
17039bd66eb0SPeter Brune @*/
1704d5def619SJonas Heinzmann PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFn *projectfunc, SNESLineSearchVINormFn *normfunc, SNESLineSearchVIDirDerivFn *dirderivfunc)
1705d71ae5a4SJacob Faibussowitsch {
17069bd66eb0SPeter Brune   PetscFunctionBegin;
1707f1c6b773SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
17089bd66eb0SPeter Brune   if (projectfunc) linesearch->ops->viproject = projectfunc;
17099bd66eb0SPeter Brune   if (normfunc) linesearch->ops->vinorm = normfunc;
1710d5def619SJonas Heinzmann   if (dirderivfunc) linesearch->ops->vidirderiv = dirderivfunc;
17113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17129bd66eb0SPeter Brune }
17139bd66eb0SPeter Brune 
17149bd66eb0SPeter Brune /*@C
1715f1c6b773SPeter Brune   SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.
17169bd66eb0SPeter Brune 
1717f6dfbefdSBarry Smith   Not Collective
1718f6dfbefdSBarry Smith 
1719f899ff85SJose E. Roman   Input Parameter:
1720f6dfbefdSBarry Smith . linesearch - the line search context, obtain with `SNESGetLineSearch()`
17219bd66eb0SPeter Brune 
17229bd66eb0SPeter Brune   Output Parameters:
17238434afd1SBarry Smith + projectfunc  - function for projecting the function to the bounds, see `SNESLineSearchVIProjectFn` for calling sequence
1724d5def619SJonas Heinzmann . normfunc     - function for computing the norm of an active set, see `SNESLineSearchVINormFn ` for calling sequence
1725d5def619SJonas Heinzmann - dirderivfunc - function for computing the directional derivative of an active set, see `SNESLineSearchVIDirDerivFn` for calling sequence
17269bd66eb0SPeter Brune 
1727f6dfbefdSBarry Smith   Level: advanced
17289bd66eb0SPeter Brune 
17299bcc50f1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetVIFunctions()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()`,
17308434afd1SBarry Smith           `SNESLineSearchVIProjectFn`, `SNESLineSearchVINormFn`
17319bd66eb0SPeter Brune @*/
1732d5def619SJonas Heinzmann PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFn **projectfunc, SNESLineSearchVINormFn **normfunc, SNESLineSearchVIDirDerivFn **dirderivfunc)
1733d71ae5a4SJacob Faibussowitsch {
17349bd66eb0SPeter Brune   PetscFunctionBegin;
17359bd66eb0SPeter Brune   if (projectfunc) *projectfunc = linesearch->ops->viproject;
17369bd66eb0SPeter Brune   if (normfunc) *normfunc = linesearch->ops->vinorm;
1737d5def619SJonas Heinzmann   if (dirderivfunc) *dirderivfunc = linesearch->ops->vidirderiv;
17383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17399bd66eb0SPeter Brune }
17409bd66eb0SPeter Brune 
1741bf7f4e0aSPeter Brune /*@C
1742420bcc1bSBarry Smith   SNESLineSearchRegister - register a line search type `SNESLineSearchType`
1743ceaaa498SBarry Smith 
1744cc4c1da9SBarry Smith   Logically Collective, No Fortran Support
1745cc4c1da9SBarry Smith 
1746ceaaa498SBarry Smith   Input Parameters:
1747420bcc1bSBarry Smith + sname    - name of the `SNESLineSearchType()`
1748ceaaa498SBarry Smith - function - the creation function for that type
1749ceaaa498SBarry Smith 
1750ceaaa498SBarry Smith   Calling sequence of `function`:
1751ceaaa498SBarry Smith . ls - the line search context
1752bf7f4e0aSPeter Brune 
1753bf7f4e0aSPeter Brune   Level: advanced
1754f6dfbefdSBarry Smith 
1755420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchSetType()`
1756bf7f4e0aSPeter Brune @*/
1757ceaaa498SBarry Smith PetscErrorCode SNESLineSearchRegister(const char sname[], PetscErrorCode (*function)(SNESLineSearch ls))
1758d71ae5a4SJacob Faibussowitsch {
1759bf7f4e0aSPeter Brune   PetscFunctionBegin;
17609566063dSJacob Faibussowitsch   PetscCall(SNESInitializePackage());
17619566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&SNESLineSearchList, sname, function));
17623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1763bf7f4e0aSPeter Brune }
1764