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 @*/
SNESLineSearchMonitorCancel(SNESLineSearch ls)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 @*/
SNESLineSearchMonitor(SNESLineSearch ls)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 @*/
SNESLineSearchMonitorSet(SNESLineSearch ls,PetscErrorCode (* f)(SNESLineSearch ls,PetscCtx mctx),PetscCtx mctx,PetscCtxDestroyFn * monitordestroy)99*2a8381b2SBarry Smith PetscErrorCode SNESLineSearchMonitorSet(SNESLineSearch ls, PetscErrorCode (*f)(SNESLineSearch ls, PetscCtx mctx), PetscCtx 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 @*/
SNESLineSearchMonitorSolutionUpdate(SNESLineSearch ls,PetscViewerAndFormat * vf)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 @*/
SNESLineSearchCreate(MPI_Comm comm,SNESLineSearch * outlinesearch)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;
19076c63389SBarry Smith linesearch->reason = 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 @*/
SNESLineSearchSetUp(SNESLineSearch linesearch)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 @*/
SNESLineSearchReset(SNESLineSearch linesearch)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 @*/
SNESLineSearchSetFunction(SNESLineSearch linesearch,PetscErrorCode (* func)(SNES snes,Vec x,Vec f))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 @*/
SNESLineSearchSetPreCheck(SNESLineSearch linesearch,PetscErrorCode (* func)(SNESLineSearch ls,Vec x,Vec d,PetscBool * changed_d,PetscCtx ctx),PetscCtx ctx)331*2a8381b2SBarry Smith PetscErrorCode SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch ls, Vec x, Vec d, PetscBool *changed_d, PetscCtx ctx), PetscCtx 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 @*/
SNESLineSearchGetPreCheck(SNESLineSearch linesearch,PetscErrorCode (** func)(SNESLineSearch,Vec,Vec,PetscBool *,void *),PetscCtxRt ctx)354*2a8381b2SBarry Smith PetscErrorCode SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch, Vec, Vec, PetscBool *, void *), PetscCtxRt ctx)
355d71ae5a4SJacob Faibussowitsch {
3569bd66eb0SPeter Brune PetscFunctionBegin;
357f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
358f190f2fcSBarry Smith if (func) *func = linesearch->ops->precheck;
359*2a8381b2SBarry Smith if (ctx) *(void **)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 @*/
SNESLineSearchSetPostCheck(SNESLineSearch linesearch,PetscErrorCode (* func)(SNESLineSearch ls,Vec x,Vec d,Vec w,PetscBool * changed_d,PetscBool * changed_w,PetscCtx ctx),PetscCtx ctx)394*2a8381b2SBarry Smith PetscErrorCode SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch ls, Vec x, Vec d, Vec w, PetscBool *changed_d, PetscBool *changed_w, PetscCtx ctx), PetscCtx 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 @*/
SNESLineSearchGetPostCheck(SNESLineSearch linesearch,PetscErrorCode (** func)(SNESLineSearch,Vec,Vec,Vec,PetscBool *,PetscBool *,void *),PetscCtxRt ctx)417*2a8381b2SBarry Smith PetscErrorCode SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), PetscCtxRt ctx)
418d71ae5a4SJacob Faibussowitsch {
4199bd66eb0SPeter Brune PetscFunctionBegin;
420f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
421f190f2fcSBarry Smith if (func) *func = linesearch->ops->postcheck;
422*2a8381b2SBarry Smith if (ctx) *(void **)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 @*/
SNESLineSearchPreCheck(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool * changed)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 @*/
SNESLineSearchPostCheck(SNESLineSearch linesearch,Vec X,Vec Y,Vec W,PetscBool * changed_Y,PetscBool * changed_W)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 @*/
SNESLineSearchPreCheckPicard(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool * changed,PetscCtx ctx)533*2a8381b2SBarry Smith PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed, PetscCtx 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
60776c63389SBarry Smith Level: advanced
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
61776c63389SBarry Smith In certain situations `SNESLineSearchApply()` may directly set a `SNESConvergedReason` in the `SNES` object so one should always check
61876c63389SBarry Smith this value immediately after the call to `SNESLineSearchApply()`.
61976c63389SBarry Smith
6201bd63e3eSSatish Balay .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchGetLambda()`, `SNESLineSearchPreCheck()`, `SNESLineSearchPostCheck()`, `SNESSolve()`, `SNESComputeFunction()`, `SNESLineSearchSetComputeNorms()`,
621db781477SPatrick Sanan `SNESLineSearchType`, `SNESLineSearchSetType()`
622f40b411bSPeter Brune @*/
SNESLineSearchApply(SNESLineSearch linesearch,Vec X,Vec F,PetscReal * fnorm,Vec Y)623d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal *fnorm, Vec Y)
624d71ae5a4SJacob Faibussowitsch {
625bf388a1fSBarry Smith PetscFunctionBegin;
626f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
627bf7f4e0aSPeter Brune PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
628bf7f4e0aSPeter Brune PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
629064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(Y, VEC_CLASSID, 5);
630bf7f4e0aSPeter Brune
63176c63389SBarry Smith linesearch->reason = SNES_LINESEARCH_SUCCEEDED;
632bf7f4e0aSPeter Brune
633bf7f4e0aSPeter Brune linesearch->vec_sol = X;
634bf7f4e0aSPeter Brune linesearch->vec_update = Y;
635bf7f4e0aSPeter Brune linesearch->vec_func = F;
636bf7f4e0aSPeter Brune
6379566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetUp(linesearch));
638bf7f4e0aSPeter Brune
639f5af7f23SKarl Rupp if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */
640bf7f4e0aSPeter Brune
641f5af7f23SKarl Rupp if (fnorm) linesearch->fnorm = *fnorm;
64248a46eb9SPierre Jolivet else PetscCall(VecNorm(F, NORM_2, &linesearch->fnorm));
643bf7f4e0aSPeter Brune
6449566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(SNESLINESEARCH_Apply, linesearch, X, F, Y));
645bf7f4e0aSPeter Brune
646dbbe0bcdSBarry Smith PetscUseTypeMethod(linesearch, apply);
647bf7f4e0aSPeter Brune
6489566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(SNESLINESEARCH_Apply, linesearch, X, F, Y));
649bf7f4e0aSPeter Brune
650f5af7f23SKarl Rupp if (fnorm) *fnorm = linesearch->fnorm;
6513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
652bf7f4e0aSPeter Brune }
653bf7f4e0aSPeter Brune
654f40b411bSPeter Brune /*@
655f1c6b773SPeter Brune SNESLineSearchDestroy - Destroys the line search instance.
656f40b411bSPeter Brune
657c3339decSBarry Smith Collective
658f40b411bSPeter Brune
6592fe279fdSBarry Smith Input Parameter:
6608e557f58SPeter Brune . linesearch - The line search context
661f40b411bSPeter Brune
66284238204SBarry Smith Level: developer
663f40b411bSPeter Brune
664420bcc1bSBarry Smith Note:
665420bcc1bSBarry Smith The line search in `SNES` is automatically called on `SNESDestroy()` so this call is rarely needed
666420bcc1bSBarry Smith
667420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchReset()`, `SNESDestroy()`
668f40b411bSPeter Brune @*/
SNESLineSearchDestroy(SNESLineSearch * linesearch)669d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchDestroy(SNESLineSearch *linesearch)
670d71ae5a4SJacob Faibussowitsch {
671bf7f4e0aSPeter Brune PetscFunctionBegin;
6723ba16761SJacob Faibussowitsch if (!*linesearch) PetscFunctionReturn(PETSC_SUCCESS);
673f4f49eeaSPierre Jolivet PetscValidHeaderSpecific(*linesearch, SNESLINESEARCH_CLASSID, 1);
674f4f49eeaSPierre Jolivet if (--((PetscObject)*linesearch)->refct > 0) {
6759371c9d4SSatish Balay *linesearch = NULL;
6763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
6779371c9d4SSatish Balay }
6789566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsViewOff((PetscObject)*linesearch));
6799566063dSJacob Faibussowitsch PetscCall(SNESLineSearchReset(*linesearch));
6803ba16761SJacob Faibussowitsch PetscTryTypeMethod(*linesearch, destroy);
681648c30bcSBarry Smith PetscCall(PetscViewerDestroy(&(*linesearch)->monitor));
682f4f49eeaSPierre Jolivet PetscCall(SNESLineSearchMonitorCancel(*linesearch));
6839566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(linesearch));
6843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
685bf7f4e0aSPeter Brune }
686bf7f4e0aSPeter Brune
687f40b411bSPeter Brune /*@
688dcf2fd19SBarry Smith SNESLineSearchSetDefaultMonitor - Turns on/off printing useful information and debugging output about the line search.
689bf7f4e0aSPeter Brune
690c3339decSBarry Smith Logically Collective
691f6dfbefdSBarry Smith
692bf7f4e0aSPeter Brune Input Parameters:
693dcf2fd19SBarry Smith + linesearch - the linesearch object
69420f4b53cSBarry Smith - viewer - an `PETSCVIEWERASCII` `PetscViewer` or `NULL` to turn off monitor
695bf7f4e0aSPeter Brune
696f6dfbefdSBarry Smith Options Database Key:
697dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - enables the monitor
698bf7f4e0aSPeter Brune
699bf7f4e0aSPeter Brune Level: intermediate
700bf7f4e0aSPeter Brune
701e4094ef1SJacob Faibussowitsch Developer Notes:
702f6dfbefdSBarry Smith This monitor is implemented differently than the other line search monitors that are set with
703f6dfbefdSBarry Smith `SNESLineSearchMonitorSet()` since it is called in many locations of the line search routines to display aspects of the
704d12e167eSBarry Smith line search that are not visible to the other monitors.
705dcf2fd19SBarry Smith
706420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `PETSCVIEWERASCII`, `SNESGetLineSearch()`, `SNESLineSearchGetDefaultMonitor()`, `PetscViewer`, `SNESLineSearchSetMonitor()`,
707f6dfbefdSBarry Smith `SNESLineSearchMonitorSetFromOptions()`
708bf7f4e0aSPeter Brune @*/
SNESLineSearchSetDefaultMonitor(SNESLineSearch linesearch,PetscViewer viewer)709d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetDefaultMonitor(SNESLineSearch linesearch, PetscViewer viewer)
710d71ae5a4SJacob Faibussowitsch {
711bf7f4e0aSPeter Brune PetscFunctionBegin;
712648c30bcSBarry Smith PetscCall(PetscViewerDestroy(&linesearch->monitor));
713dcf2fd19SBarry Smith linesearch->monitor = viewer;
7143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
715bf7f4e0aSPeter Brune }
716bf7f4e0aSPeter Brune
717f40b411bSPeter Brune /*@
718420bcc1bSBarry Smith SNESLineSearchGetDefaultMonitor - Gets the `PetscViewer` instance for the default line search monitor that is turned on with `SNESLineSearchSetDefaultMonitor()`
719f6dfbefdSBarry Smith
720c3339decSBarry Smith Logically Collective
7216a388c36SPeter Brune
722f190f2fcSBarry Smith Input Parameter:
723420bcc1bSBarry Smith . linesearch - the line search context
724f40b411bSPeter Brune
725f190f2fcSBarry Smith Output Parameter:
7268e557f58SPeter Brune . monitor - monitor context
727f40b411bSPeter Brune
728f40b411bSPeter Brune Level: intermediate
729f40b411bSPeter Brune
730420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchSetDefaultMonitor()`, `PetscViewer`
731f40b411bSPeter Brune @*/
SNESLineSearchGetDefaultMonitor(SNESLineSearch linesearch,PetscViewer * monitor)732d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetDefaultMonitor(SNESLineSearch linesearch, PetscViewer *monitor)
733d71ae5a4SJacob Faibussowitsch {
7346a388c36SPeter Brune PetscFunctionBegin;
735f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
7366a388c36SPeter Brune *monitor = linesearch->monitor;
7373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
7386a388c36SPeter Brune }
7396a388c36SPeter Brune
740dcf2fd19SBarry Smith /*@C
741420bcc1bSBarry Smith SNESLineSearchMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated in the options database
742dcf2fd19SBarry Smith
743c3339decSBarry Smith Collective
744dcf2fd19SBarry Smith
745dcf2fd19SBarry Smith Input Parameters:
746420bcc1bSBarry Smith + ls - `SNESLineSearch` object to monitor
747420bcc1bSBarry Smith . name - the monitor type
748dcf2fd19SBarry Smith . help - message indicating what monitoring is done
749dcf2fd19SBarry Smith . manual - manual page for the monitor
75049abdd8aSBarry Smith . monitor - the monitor function, must use `PetscViewerAndFormat` as its context
751f6dfbefdSBarry 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`
752dcf2fd19SBarry Smith
753420bcc1bSBarry Smith Calling sequence of `monitor`:
754420bcc1bSBarry Smith + ls - `SNESLineSearch` object being monitored
755420bcc1bSBarry Smith - vf - a `PetscViewerAndFormat` struct that provides the `PetscViewer` and `PetscViewerFormat` being used
756dcf2fd19SBarry Smith
757420bcc1bSBarry Smith Calling sequence of `monitorsetup`:
758420bcc1bSBarry Smith + ls - `SNESLineSearch` object being monitored
759420bcc1bSBarry Smith - vf - a `PetscViewerAndFormat` struct that provides the `PetscViewer` and `PetscViewerFormat` being used
760420bcc1bSBarry Smith
761420bcc1bSBarry Smith Level: advanced
762420bcc1bSBarry Smith
763648c30bcSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetMonitor()`, `PetscOptionsCreateViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
764db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
765e4094ef1SJacob Faibussowitsch `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
766db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
767c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
768db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
769db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()`
770dcf2fd19SBarry Smith @*/
SNESLineSearchMonitorSetFromOptions(SNESLineSearch ls,const char name[],const char help[],const char manual[],PetscErrorCode (* monitor)(SNESLineSearch ls,PetscViewerAndFormat * vf),PetscErrorCode (* monitorsetup)(SNESLineSearch ls,PetscViewerAndFormat * vf))771420bcc1bSBarry 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))
772d71ae5a4SJacob Faibussowitsch {
773dcf2fd19SBarry Smith PetscViewer viewer;
774dcf2fd19SBarry Smith PetscViewerFormat format;
775dcf2fd19SBarry Smith PetscBool flg;
776dcf2fd19SBarry Smith
777dcf2fd19SBarry Smith PetscFunctionBegin;
778648c30bcSBarry Smith PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)ls), ((PetscObject)ls)->options, ((PetscObject)ls)->prefix, name, &viewer, &format, &flg));
779dcf2fd19SBarry Smith if (flg) {
780d12e167eSBarry Smith PetscViewerAndFormat *vf;
7819566063dSJacob Faibussowitsch PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf));
782648c30bcSBarry Smith PetscCall(PetscViewerDestroy(&viewer));
7831baa6e33SBarry Smith if (monitorsetup) PetscCall((*monitorsetup)(ls, vf));
784*2a8381b2SBarry Smith PetscCall(SNESLineSearchMonitorSet(ls, (PetscErrorCode (*)(SNESLineSearch, PetscCtx))monitor, vf, (PetscCtxDestroyFn *)PetscViewerAndFormatDestroy));
785dcf2fd19SBarry Smith }
7863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
787dcf2fd19SBarry Smith }
788dcf2fd19SBarry Smith
789f40b411bSPeter Brune /*@
790f1c6b773SPeter Brune SNESLineSearchSetFromOptions - Sets options for the line search
791f40b411bSPeter Brune
792c3339decSBarry Smith Logically Collective
793f6dfbefdSBarry Smith
794f899ff85SJose E. Roman Input Parameter:
795420bcc1bSBarry Smith . linesearch - a `SNESLineSearch` line search context
796f40b411bSPeter Brune
797f40b411bSPeter Brune Options Database Keys:
798a99ef635SJonas Heinzmann + -snes_linesearch_type <type> - basic (or equivalently none), `bt`, `secant`, `cp`, `nleqerr`, `bisection`, `shell`
799a99ef635SJonas Heinzmann . -snes_linesearch_order <order> - 1, 2, 3. Most types only support certain orders (`bt` supports 1, 2 or 3)
800f6dfbefdSBarry Smith . -snes_linesearch_norms - Turn on/off the linesearch norms for the basic linesearch typem (`SNESLineSearchSetComputeNorms()`)
801a99ef635SJonas Heinzmann . -snes_linesearch_minlambda - The minimum `lambda`
802a99ef635SJonas Heinzmann . -snes_linesearch_maxlambda - The maximum `lambda`
8031a9b3a06SPeter Brune . -snes_linesearch_rtol - Relative tolerance for iterative line searches
8041a9b3a06SPeter Brune . -snes_linesearch_atol - Absolute tolerance for iterative line searches
805a99ef635SJonas Heinzmann . -snes_linesearch_ltol - Change in `lambda` tolerance for iterative line searches
8061a9b3a06SPeter Brune . -snes_linesearch_max_it - The number of iterations for iterative line searches
807dcf2fd19SBarry Smith . -snes_linesearch_monitor [:filename] - Print progress of line searches
808dcf2fd19SBarry Smith . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine
8098e557f58SPeter Brune . -snes_linesearch_damping - The linesearch damping parameter
810a99ef635SJonas Heinzmann . -snes_linesearch_keeplambda - Keep the previous `lambda` as the initial guess.
8111a9b3a06SPeter Brune . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method
812d8d34be6SBarry Smith - -snes_linesearch_precheck_picard_angle - Angle used in Picard precheck method
813f40b411bSPeter Brune
814f40b411bSPeter Brune Level: intermediate
815f40b411bSPeter Brune
81662842358SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchSetOrder()`, `SNESLineSearchSetType()`, `SNESLineSearchSetTolerances()`, `SNESLineSearchSetDamping()`, `SNESLineSearchPreCheckPicard()`,
817db781477SPatrick Sanan `SNESLineSearchType`, `SNESLineSearchSetComputeNorms()`
818f40b411bSPeter Brune @*/
SNESLineSearchSetFromOptions(SNESLineSearch linesearch)819d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch)
820d71ae5a4SJacob Faibussowitsch {
8211a4f838cSPeter Brune const char *deft = SNESLINESEARCHBASIC;
822bf7f4e0aSPeter Brune char type[256];
823bf7f4e0aSPeter Brune PetscBool flg, set;
824dcf2fd19SBarry Smith PetscViewer viewer;
825bf388a1fSBarry Smith
826bf7f4e0aSPeter Brune PetscFunctionBegin;
8279566063dSJacob Faibussowitsch PetscCall(SNESLineSearchRegisterAll());
828bf7f4e0aSPeter Brune
829d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)linesearch);
830f5af7f23SKarl Rupp if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name;
8319566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-snes_linesearch_type", "Linesearch type", "SNESLineSearchSetType", SNESLineSearchList, deft, type, 256, &flg));
832bf7f4e0aSPeter Brune if (flg) {
8339566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch, type));
834bf7f4e0aSPeter Brune } else if (!((PetscObject)linesearch)->type_name) {
8359566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch, deft));
836bf7f4e0aSPeter Brune }
837bf7f4e0aSPeter Brune
838648c30bcSBarry Smith PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)linesearch), ((PetscObject)linesearch)->options, ((PetscObject)linesearch)->prefix, "-snes_linesearch_monitor", &viewer, NULL, &set));
839648c30bcSBarry Smith if (set) PetscCall(SNESLineSearchSetDefaultMonitor(linesearch, viewer));
8409566063dSJacob Faibussowitsch PetscCall(SNESLineSearchMonitorSetFromOptions(linesearch, "-snes_linesearch_monitor_solution_update", "View correction at each iteration", "SNESLineSearchMonitorSolutionUpdate", SNESLineSearchMonitorSolutionUpdate, NULL));
841bf7f4e0aSPeter Brune
8421a9b3a06SPeter Brune /* tolerances */
843a99ef635SJonas Heinzmann PetscCall(PetscOptionsReal("-snes_linesearch_minlambda", "Minimum lambda", "SNESLineSearchSetTolerances", linesearch->minlambda, &linesearch->minlambda, NULL));
844a99ef635SJonas Heinzmann PetscCall(PetscOptionsReal("-snes_linesearch_maxlambda", "Maximum lambda", "SNESLineSearchSetTolerances", linesearch->maxlambda, &linesearch->maxlambda, NULL));
8459566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_linesearch_rtol", "Relative tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->rtol, &linesearch->rtol, NULL));
8469566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_linesearch_atol", "Absolute tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->atol, &linesearch->atol, NULL));
8479566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_linesearch_ltol", "Change in lambda tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->ltol, &linesearch->ltol, NULL));
848a99ef635SJonas Heinzmann PetscCall(PetscOptionsInt("-snes_linesearch_max_it", "Maximum iterations for iterative line searches", "SNESLineSearchSetTolerances", linesearch->max_it, &linesearch->max_it, NULL));
849a99ef635SJonas Heinzmann
850a99ef635SJonas Heinzmann /* deprecated options */
851a99ef635SJonas Heinzmann PetscCall(PetscOptionsDeprecated("-snes_linesearch_maxstep", "-snes_linesearch_maxlambda", "3.24.0", NULL));
8527a35526eSPeter Brune
8531a9b3a06SPeter Brune /* damping parameters */
854a99ef635SJonas Heinzmann PetscCall(PetscOptionsReal("-snes_linesearch_damping", "Line search damping (and depending on chosen line search initial lambda guess)", "SNESLineSearchSetDamping", linesearch->damping, &linesearch->damping, NULL));
8551a9b3a06SPeter Brune
8569566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_linesearch_keeplambda", "Use previous lambda as damping", "SNESLineSearchSetKeepLambda", linesearch->keeplambda, &linesearch->keeplambda, NULL));
8571a9b3a06SPeter Brune
8581a9b3a06SPeter Brune /* precheck */
8599566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_linesearch_precheck_picard", "Use a correction that sometimes improves convergence of Picard iteration", "SNESLineSearchPreCheckPicard", flg, &flg, &set));
8607a35526eSPeter Brune if (set) {
8617a35526eSPeter Brune if (flg) {
8627a35526eSPeter Brune linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */
863f5af7f23SKarl Rupp
864d0609cedSBarry 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));
8659566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetPreCheck(linesearch, SNESLineSearchPreCheckPicard, &linesearch->precheck_picard_angle));
8667a35526eSPeter Brune } else {
8679566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetPreCheck(linesearch, NULL, NULL));
8687a35526eSPeter Brune }
8697a35526eSPeter Brune }
8709566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_linesearch_order", "Order of approximation used in the line search", "SNESLineSearchSetOrder", linesearch->order, &linesearch->order, NULL));
8719566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_linesearch_norms", "Compute final norms in line search", "SNESLineSearchSetComputeNorms", linesearch->norms, &linesearch->norms, NULL));
8727a35526eSPeter Brune
873dbbe0bcdSBarry Smith PetscTryTypeMethod(linesearch, setfromoptions, PetscOptionsObject);
8745a9b6599SPeter Brune
875dbbe0bcdSBarry Smith PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)linesearch, PetscOptionsObject));
876d0609cedSBarry Smith PetscOptionsEnd();
8773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
878bf7f4e0aSPeter Brune }
879bf7f4e0aSPeter Brune
880f40b411bSPeter Brune /*@
881f190f2fcSBarry Smith SNESLineSearchView - Prints useful information about the line search
882f40b411bSPeter Brune
88320f4b53cSBarry Smith Logically Collective
88420f4b53cSBarry Smith
885f40b411bSPeter Brune Input Parameters:
8862fe279fdSBarry Smith + linesearch - line search context
887420bcc1bSBarry Smith - viewer - the `PetscViewer` to display the line search information to
888f40b411bSPeter Brune
889f40b411bSPeter Brune Level: intermediate
890f40b411bSPeter Brune
891420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `PetscViewer`, `SNESLineSearchCreate()`
892f40b411bSPeter Brune @*/
SNESLineSearchView(SNESLineSearch linesearch,PetscViewer viewer)893d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer)
894d71ae5a4SJacob Faibussowitsch {
8959f196a02SMartin Diehl PetscBool isascii;
896bf388a1fSBarry Smith
897bf7f4e0aSPeter Brune PetscFunctionBegin;
8987f1410a3SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
89948a46eb9SPierre Jolivet if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch), &viewer));
9007f1410a3SPeter Brune PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
9017f1410a3SPeter Brune PetscCheckSameComm(linesearch, 1, viewer, 2);
902f40b411bSPeter Brune
9039f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
9049f196a02SMartin Diehl if (isascii) {
9059566063dSJacob Faibussowitsch PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)linesearch, viewer));
9069566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer));
907dbbe0bcdSBarry Smith PetscTryTypeMethod(linesearch, view, viewer);
9089566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer));
909a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIIPrintf(viewer, " maxlambda=%e, minlambda=%e\n", (double)linesearch->maxlambda, (double)linesearch->minlambda));
9109566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol, (double)linesearch->atol, (double)linesearch->ltol));
911a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIIPrintf(viewer, " maximum iterations=%" PetscInt_FMT "\n", linesearch->max_it));
9126b2b7091SBarry Smith if (linesearch->ops->precheck) {
9136b2b7091SBarry Smith if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) {
91463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " using precheck step to speed up Picard convergence\n"));
9157f1410a3SPeter Brune } else {
91663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " using user-defined precheck step\n"));
9177f1410a3SPeter Brune }
9187f1410a3SPeter Brune }
91948a46eb9SPierre Jolivet if (linesearch->ops->postcheck) PetscCall(PetscViewerASCIIPrintf(viewer, " using user-defined postcheck step\n"));
9207f1410a3SPeter Brune }
9213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
922bf7f4e0aSPeter Brune }
923bf7f4e0aSPeter Brune
924cc4c1da9SBarry Smith /*@
925420bcc1bSBarry Smith SNESLineSearchGetType - Gets the `SNESLinesearchType` of a `SNESLineSearch`
926a80ff896SJed Brown
927c3339decSBarry Smith Logically Collective
928a80ff896SJed Brown
9292fe279fdSBarry Smith Input Parameter:
930420bcc1bSBarry Smith . linesearch - the line search context
931a80ff896SJed Brown
9322fe279fdSBarry Smith Output Parameter:
9332fe279fdSBarry Smith . type - The type of line search, or `NULL` if not set
934a80ff896SJed Brown
935a80ff896SJed Brown Level: intermediate
936a80ff896SJed Brown
937420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchSetType()`
938a80ff896SJed Brown @*/
SNESLineSearchGetType(SNESLineSearch linesearch,SNESLineSearchType * type)939d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetType(SNESLineSearch linesearch, SNESLineSearchType *type)
940d71ae5a4SJacob Faibussowitsch {
941a80ff896SJed Brown PetscFunctionBegin;
942a80ff896SJed Brown PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
9434f572ea9SToby Isaac PetscAssertPointer(type, 2);
944a80ff896SJed Brown *type = ((PetscObject)linesearch)->type_name;
9453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
946a80ff896SJed Brown }
947a80ff896SJed Brown
948cc4c1da9SBarry Smith /*@
9490b4b7b1cSBarry Smith SNESLineSearchSetType - Sets the `SNESLinesearchType` of a `SNESLineSearch` object to indicate the line search algorithm that should be used by a given `SNES` solver
950f40b411bSPeter Brune
951c3339decSBarry Smith Logically Collective
952f190f2fcSBarry Smith
953f40b411bSPeter Brune Input Parameters:
954420bcc1bSBarry Smith + linesearch - the line search context
955ceaaa498SBarry Smith - type - The type of line search to be used, see `SNESLineSearchType`
9561fe24845SBarry Smith
9573c7db156SBarry Smith Options Database Key:
958a99ef635SJonas Heinzmann . -snes_linesearch_type <type> - basic (or equivalently none), bt, secant, cp, nleqerr, bisection, shell
959cd7522eaSPeter Brune
960f40b411bSPeter Brune Level: intermediate
961f40b411bSPeter Brune
9620b4b7b1cSBarry Smith Note:
9630b4b7b1cSBarry Smith The `SNESLineSearch` object is generally obtained with `SNESGetLineSearch()`
9640b4b7b1cSBarry Smith
9650b4b7b1cSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchGetType()`,
9660b4b7b1cSBarry Smith `SNESGetLineSearch()`
967f40b411bSPeter Brune @*/
SNESLineSearchSetType(SNESLineSearch linesearch,SNESLineSearchType type)968d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type)
969d71ae5a4SJacob Faibussowitsch {
970bf7f4e0aSPeter Brune PetscBool match;
9715f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(SNESLineSearch);
972bf7f4e0aSPeter Brune
973bf7f4e0aSPeter Brune PetscFunctionBegin;
974f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
9754f572ea9SToby Isaac PetscAssertPointer(type, 2);
976bf7f4e0aSPeter Brune
9779566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)linesearch, type, &match));
9783ba16761SJacob Faibussowitsch if (match) PetscFunctionReturn(PETSC_SUCCESS);
979bf7f4e0aSPeter Brune
9809566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(SNESLineSearchList, type, &r));
9816adde796SStefano Zampini PetscCheck(r, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested Line Search type %s", type);
982bf7f4e0aSPeter Brune /* Destroy the previous private line search context */
983213acdd3SPierre Jolivet PetscTryTypeMethod(linesearch, destroy);
9840298fd71SBarry Smith linesearch->ops->destroy = NULL;
985f1c6b773SPeter Brune /* Reinitialize function pointers in SNESLineSearchOps structure */
9869e5d0892SLisandro Dalcin linesearch->ops->apply = NULL;
9879e5d0892SLisandro Dalcin linesearch->ops->view = NULL;
9889e5d0892SLisandro Dalcin linesearch->ops->setfromoptions = NULL;
9899e5d0892SLisandro Dalcin linesearch->ops->destroy = NULL;
990bf7f4e0aSPeter Brune
9919566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)linesearch, type));
9929566063dSJacob Faibussowitsch PetscCall((*r)(linesearch));
9933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
994bf7f4e0aSPeter Brune }
995bf7f4e0aSPeter Brune
996f40b411bSPeter Brune /*@
997f6dfbefdSBarry Smith SNESLineSearchSetSNES - Sets the `SNES` for the linesearch for function evaluation.
998f40b411bSPeter Brune
999f40b411bSPeter Brune Input Parameters:
1000420bcc1bSBarry Smith + linesearch - the line search context
1001420bcc1bSBarry Smith - snes - The `SNES` instance
1002f40b411bSPeter Brune
100378bcb3b5SPeter Brune Level: developer
100478bcb3b5SPeter Brune
1005f6dfbefdSBarry Smith Note:
1006f190f2fcSBarry Smith This happens automatically when the line search is obtained/created with
1007f6dfbefdSBarry Smith `SNESGetLineSearch()`. This routine is therefore mainly called within `SNES`
100878bcb3b5SPeter Brune implementations.
1009f40b411bSPeter Brune
1010420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetSNES()`, `SNESLineSearchSetVecs()`
1011f40b411bSPeter Brune @*/
SNESLineSearchSetSNES(SNESLineSearch linesearch,SNES snes)1012d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes)
1013d71ae5a4SJacob Faibussowitsch {
1014bf7f4e0aSPeter Brune PetscFunctionBegin;
1015f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1016bf7f4e0aSPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 2);
1017bf7f4e0aSPeter Brune linesearch->snes = snes;
10183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1019bf7f4e0aSPeter Brune }
1020bf7f4e0aSPeter Brune
1021f40b411bSPeter Brune /*@
1022f6dfbefdSBarry Smith SNESLineSearchGetSNES - Gets the `SNES` instance associated with the line search.
1023f6dfbefdSBarry Smith
1024f6dfbefdSBarry Smith Not Collective
1025f40b411bSPeter Brune
10262fe279fdSBarry Smith Input Parameter:
1027420bcc1bSBarry Smith . linesearch - the line search context
1028f40b411bSPeter Brune
10292fe279fdSBarry Smith Output Parameter:
10302fe279fdSBarry Smith . snes - The `SNES` instance
1031f40b411bSPeter Brune
10328141a3b9SPeter Brune Level: developer
1033f40b411bSPeter Brune
1034420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESType`, `SNESLineSearchSetVecs()`
1035f40b411bSPeter Brune @*/
SNESLineSearchGetSNES(SNESLineSearch linesearch,SNES * snes)1036d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes)
1037d71ae5a4SJacob Faibussowitsch {
1038bf7f4e0aSPeter Brune PetscFunctionBegin;
1039f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
10404f572ea9SToby Isaac PetscAssertPointer(snes, 2);
1041bf7f4e0aSPeter Brune *snes = linesearch->snes;
10423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1043bf7f4e0aSPeter Brune }
1044bf7f4e0aSPeter Brune
1045f40b411bSPeter Brune /*@
1046a99ef635SJonas Heinzmann SNESLineSearchGetLambda - Gets the last line search `lambda` used
1047f40b411bSPeter Brune
1048f6dfbefdSBarry Smith Not Collective
1049f6dfbefdSBarry Smith
10502fe279fdSBarry Smith Input Parameter:
1051420bcc1bSBarry Smith . linesearch - the line search context
1052f40b411bSPeter Brune
10532fe279fdSBarry Smith Output Parameter:
10548c5add6aSPierre Jolivet . lambda - The last `lambda` (scaling of the solution update) computed during `SNESLineSearchApply()`
1055f40b411bSPeter Brune
105678bcb3b5SPeter Brune Level: advanced
105778bcb3b5SPeter Brune
1058f6dfbefdSBarry Smith Note:
10598e557f58SPeter Brune This is useful in methods where the solver is ill-scaled and
106078bcb3b5SPeter Brune requires some adaptive notion of the difference in scale between the
1061f6dfbefdSBarry Smith solution and the function. For instance, `SNESQN` may be scaled by the
1062a99ef635SJonas Heinzmann line search `lambda` using the argument -snes_qn_scaling ls.
106378bcb3b5SPeter Brune
1064420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetLambda()`, `SNESLineSearchGetDamping()`, `SNESLineSearchApply()`
1065f40b411bSPeter Brune @*/
SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal * lambda)1066d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch linesearch, PetscReal *lambda)
1067d71ae5a4SJacob Faibussowitsch {
10686a388c36SPeter Brune PetscFunctionBegin;
1069f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
10704f572ea9SToby Isaac PetscAssertPointer(lambda, 2);
10716a388c36SPeter Brune *lambda = linesearch->lambda;
10723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
10736a388c36SPeter Brune }
10746a388c36SPeter Brune
1075f40b411bSPeter Brune /*@
1076a99ef635SJonas Heinzmann SNESLineSearchSetLambda - Sets the line search `lambda` (scaling of the solution update)
1077f40b411bSPeter Brune
1078f40b411bSPeter Brune Input Parameters:
10798e557f58SPeter Brune + linesearch - line search context
1080a99ef635SJonas Heinzmann - lambda - The `lambda` to use
1081f40b411bSPeter Brune
108220f4b53cSBarry Smith Level: advanced
108320f4b53cSBarry Smith
1084f6dfbefdSBarry Smith Note:
1085f6dfbefdSBarry Smith This routine is typically used within implementations of `SNESLineSearchApply()`
1086a99ef635SJonas Heinzmann to set the final `lambda`. This routine (and `SNESLineSearchGetLambda()`) were
1087a99ef635SJonas Heinzmann added to facilitate Quasi-Newton methods that use the previous `lambda`
1088cd7522eaSPeter Brune as an inner scaling parameter.
1089cd7522eaSPeter Brune
1090420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetLambda()`
1091f40b411bSPeter Brune @*/
SNESLineSearchSetLambda(SNESLineSearch linesearch,PetscReal lambda)1092d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
1093d71ae5a4SJacob Faibussowitsch {
10946a388c36SPeter Brune PetscFunctionBegin;
1095f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
10966a388c36SPeter Brune linesearch->lambda = lambda;
10973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
10986a388c36SPeter Brune }
10996a388c36SPeter Brune
1100f40b411bSPeter Brune /*@
1101ceaaa498SBarry Smith SNESLineSearchGetTolerances - Gets the tolerances for the line search.
1102f40b411bSPeter Brune
1103f6dfbefdSBarry Smith Not Collective
1104f6dfbefdSBarry Smith
1105f899ff85SJose E. Roman Input Parameter:
1106420bcc1bSBarry Smith . linesearch - the line search context
1107f40b411bSPeter Brune
1108f40b411bSPeter Brune Output Parameters:
1109a99ef635SJonas Heinzmann + minlambda - The minimum `lambda` allowed
1110a99ef635SJonas Heinzmann . maxlambda - The maximum `lambda` allowed
1111516fe3c3SPeter Brune . rtol - The relative tolerance for iterative line searches
1112516fe3c3SPeter Brune . atol - The absolute tolerance for iterative line searches
1113a99ef635SJonas Heinzmann . ltol - The change in `lambda` tolerance for iterative line searches
1114a99ef635SJonas Heinzmann - max_it - The maximum number of iterations of the line search
1115f40b411bSPeter Brune
111678bcb3b5SPeter Brune Level: intermediate
111778bcb3b5SPeter Brune
1118f6dfbefdSBarry Smith Note:
111978bcb3b5SPeter Brune Different line searches may implement these parameters slightly differently as
11203c7d6663SPeter Brune the type requires.
1121516fe3c3SPeter Brune
1122420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetTolerances()`
1123f40b411bSPeter Brune @*/
SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal * minlambda,PetscReal * maxlambda,PetscReal * rtol,PetscReal * atol,PetscReal * ltol,PetscInt * max_it)1124a99ef635SJonas Heinzmann PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch linesearch, PetscReal *minlambda, PetscReal *maxlambda, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_it)
1125d71ae5a4SJacob Faibussowitsch {
11266a388c36SPeter Brune PetscFunctionBegin;
1127f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1128a99ef635SJonas Heinzmann if (minlambda) {
1129a99ef635SJonas Heinzmann PetscAssertPointer(minlambda, 2);
1130a99ef635SJonas Heinzmann *minlambda = linesearch->minlambda;
1131516fe3c3SPeter Brune }
1132a99ef635SJonas Heinzmann if (maxlambda) {
1133a99ef635SJonas Heinzmann PetscAssertPointer(maxlambda, 3);
1134a99ef635SJonas Heinzmann *maxlambda = linesearch->maxlambda;
1135516fe3c3SPeter Brune }
1136516fe3c3SPeter Brune if (rtol) {
11374f572ea9SToby Isaac PetscAssertPointer(rtol, 4);
1138516fe3c3SPeter Brune *rtol = linesearch->rtol;
1139516fe3c3SPeter Brune }
1140516fe3c3SPeter Brune if (atol) {
11414f572ea9SToby Isaac PetscAssertPointer(atol, 5);
1142516fe3c3SPeter Brune *atol = linesearch->atol;
1143516fe3c3SPeter Brune }
1144516fe3c3SPeter Brune if (ltol) {
11454f572ea9SToby Isaac PetscAssertPointer(ltol, 6);
1146516fe3c3SPeter Brune *ltol = linesearch->ltol;
1147516fe3c3SPeter Brune }
1148a99ef635SJonas Heinzmann if (max_it) {
1149a99ef635SJonas Heinzmann PetscAssertPointer(max_it, 7);
1150a99ef635SJonas Heinzmann *max_it = linesearch->max_it;
1151516fe3c3SPeter Brune }
11523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
11536a388c36SPeter Brune }
11546a388c36SPeter Brune
1155f40b411bSPeter Brune /*@
1156ceaaa498SBarry Smith SNESLineSearchSetTolerances - Sets the tolerances for the linesearch.
1157f40b411bSPeter Brune
1158f6dfbefdSBarry Smith Collective
1159f6dfbefdSBarry Smith
1160f40b411bSPeter Brune Input Parameters:
1161420bcc1bSBarry Smith + linesearch - the line search context
1162a99ef635SJonas Heinzmann . minlambda - The minimum `lambda` allowed
1163a99ef635SJonas Heinzmann . maxlambda - The maximum `lambda` allowed
1164516fe3c3SPeter Brune . rtol - The relative tolerance for iterative line searches
1165516fe3c3SPeter Brune . atol - The absolute tolerance for iterative line searches
1166a99ef635SJonas Heinzmann . ltol - The change in `lambda` tolerance for iterative line searches
1167420bcc1bSBarry Smith - max_it - The maximum number of iterations of the line search
1168420bcc1bSBarry Smith
1169420bcc1bSBarry Smith Options Database Keys:
1170a99ef635SJonas Heinzmann + -snes_linesearch_minlambda - The minimum `lambda` allowed
1171a99ef635SJonas Heinzmann . -snes_linesearch_maxlambda - The maximum `lambda` allowed
1172420bcc1bSBarry Smith . -snes_linesearch_rtol - Relative tolerance for iterative line searches
1173420bcc1bSBarry Smith . -snes_linesearch_atol - Absolute tolerance for iterative line searches
1174a99ef635SJonas Heinzmann . -snes_linesearch_ltol - Change in `lambda` tolerance for iterative line searches
1175420bcc1bSBarry Smith - -snes_linesearch_max_it - The number of iterations for iterative line searches
1176f40b411bSPeter Brune
117720f4b53cSBarry Smith Level: intermediate
117820f4b53cSBarry Smith
1179f6dfbefdSBarry Smith Note:
1180420bcc1bSBarry Smith The user may choose to not set any of the tolerances using `PETSC_DEFAULT` in place of an argument.
1181f40b411bSPeter Brune
1182420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetTolerances()`
1183f40b411bSPeter Brune @*/
SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal minlambda,PetscReal maxlambda,PetscReal rtol,PetscReal atol,PetscReal ltol,PetscInt max_it)1184a99ef635SJonas Heinzmann PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch linesearch, PetscReal minlambda, PetscReal maxlambda, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_it)
1185d71ae5a4SJacob Faibussowitsch {
11866a388c36SPeter Brune PetscFunctionBegin;
1187f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1188a99ef635SJonas Heinzmann PetscValidLogicalCollectiveReal(linesearch, minlambda, 2);
1189a99ef635SJonas Heinzmann PetscValidLogicalCollectiveReal(linesearch, maxlambda, 3);
1190d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch, rtol, 4);
1191d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch, atol, 5);
1192d3952184SSatish Balay PetscValidLogicalCollectiveReal(linesearch, ltol, 6);
1193420bcc1bSBarry Smith PetscValidLogicalCollectiveInt(linesearch, max_it, 7);
1194d3952184SSatish Balay
1195a99ef635SJonas Heinzmann if (minlambda != (PetscReal)PETSC_DEFAULT) {
1196a99ef635SJonas Heinzmann PetscCheck(minlambda >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Minimum lambda %14.12e must be non-negative", (double)minlambda);
1197a99ef635SJonas 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);
1198a99ef635SJonas Heinzmann linesearch->minlambda = minlambda;
1199d3952184SSatish Balay }
1200d3952184SSatish Balay
1201a99ef635SJonas Heinzmann if (maxlambda != (PetscReal)PETSC_DEFAULT) {
1202a99ef635SJonas Heinzmann PetscCheck(maxlambda > 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum lambda %14.12e must be positive", (double)maxlambda);
1203a99ef635SJonas Heinzmann linesearch->maxlambda = maxlambda;
1204d3952184SSatish Balay }
1205d3952184SSatish Balay
120613bcc0bdSJacob Faibussowitsch if (rtol != (PetscReal)PETSC_DEFAULT) {
12072061ca32SJunchao 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);
1208516fe3c3SPeter Brune linesearch->rtol = rtol;
1209d3952184SSatish Balay }
1210d3952184SSatish Balay
121113bcc0bdSJacob Faibussowitsch if (atol != (PetscReal)PETSC_DEFAULT) {
12125f80ce2aSJacob Faibussowitsch PetscCheck(atol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Absolute tolerance %14.12e must be non-negative", (double)atol);
1213516fe3c3SPeter Brune linesearch->atol = atol;
1214d3952184SSatish Balay }
1215d3952184SSatish Balay
121613bcc0bdSJacob Faibussowitsch if (ltol != (PetscReal)PETSC_DEFAULT) {
12175f80ce2aSJacob Faibussowitsch PetscCheck(ltol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Lambda tolerance %14.12e must be non-negative", (double)ltol);
1218516fe3c3SPeter Brune linesearch->ltol = ltol;
1219d3952184SSatish Balay }
1220d3952184SSatish Balay
1221420bcc1bSBarry Smith if (max_it != PETSC_DEFAULT) {
1222420bcc1bSBarry Smith PetscCheck(max_it >= 0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of iterations %" PetscInt_FMT " must be non-negative", max_it);
1223a99ef635SJonas Heinzmann linesearch->max_it = max_it;
1224d3952184SSatish Balay }
12253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
12266a388c36SPeter Brune }
12276a388c36SPeter Brune
1228f40b411bSPeter Brune /*@
1229f1c6b773SPeter Brune SNESLineSearchGetDamping - Gets the line search damping parameter.
1230f40b411bSPeter Brune
12312fe279fdSBarry Smith Input Parameter:
1232420bcc1bSBarry Smith . linesearch - the line search context
1233f40b411bSPeter Brune
12342fe279fdSBarry Smith Output Parameter:
12358e557f58SPeter Brune . damping - The damping parameter
1236f40b411bSPeter Brune
123778bcb3b5SPeter Brune Level: advanced
1238f40b411bSPeter Brune
1239420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearchGetStepTolerance()`, `SNESQN`
1240f40b411bSPeter Brune @*/
SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal * damping)1241d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch linesearch, PetscReal *damping)
1242d71ae5a4SJacob Faibussowitsch {
12436a388c36SPeter Brune PetscFunctionBegin;
1244f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
12454f572ea9SToby Isaac PetscAssertPointer(damping, 2);
12466a388c36SPeter Brune *damping = linesearch->damping;
12473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
12486a388c36SPeter Brune }
12496a388c36SPeter Brune
1250f40b411bSPeter Brune /*@
1251fd292e60Sprj- SNESLineSearchSetDamping - Sets the line search damping parameter.
1252f40b411bSPeter Brune
1253f40b411bSPeter Brune Input Parameters:
1254420bcc1bSBarry Smith + linesearch - the line search context
125503fd4120SBarry Smith - damping - The damping parameter
1256f40b411bSPeter Brune
1257f6dfbefdSBarry Smith Options Database Key:
1258f6dfbefdSBarry Smith . -snes_linesearch_damping <damping> - the damping value
1259f6dfbefdSBarry Smith
1260f40b411bSPeter Brune Level: intermediate
1261f40b411bSPeter Brune
1262f6dfbefdSBarry Smith Note:
1263f6dfbefdSBarry Smith The `SNESLINESEARCHNONE` line search merely takes the update step scaled by the damping parameter.
1264a99ef635SJonas Heinzmann The use of the damping parameter in the `SNESLINESEARCHSECANT` and `SNESLINESEARCHCP` line searches is much more subtle;
1265a99ef635SJonas Heinzmann it is used as a starting point for the secant method. Depending on the choice for `maxlambda`,
1266a99ef635SJonas Heinzmann the eventual `lambda` may be greater than the damping parameter however.
1267a99ef635SJonas Heinzmann For `SNESLINESEARCHBISECTION` and `SNESLINESEARCHBT` the damping is instead used as the initial guess,
1268a99ef635SJonas Heinzmann below which the line search will not go. Hence, it is the maximum possible value for `lambda`.
1269cd7522eaSPeter Brune
1270420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetDamping()`
1271f40b411bSPeter Brune @*/
SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping)1272d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch linesearch, PetscReal damping)
1273d71ae5a4SJacob Faibussowitsch {
12746a388c36SPeter Brune PetscFunctionBegin;
1275f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
12766a388c36SPeter Brune linesearch->damping = damping;
12773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
12786a388c36SPeter Brune }
12796a388c36SPeter Brune
128059405d5eSPeter Brune /*@
128159405d5eSPeter Brune SNESLineSearchGetOrder - Gets the line search approximation order.
128259405d5eSPeter Brune
1283f6dfbefdSBarry Smith Input Parameter:
1284420bcc1bSBarry Smith . linesearch - the line search context
128559405d5eSPeter Brune
1286f6dfbefdSBarry Smith Output Parameter:
12878e557f58SPeter Brune . order - The order
128859405d5eSPeter Brune
128959405d5eSPeter Brune Level: intermediate
129059405d5eSPeter Brune
1291420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetOrder()`
129259405d5eSPeter Brune @*/
SNESLineSearchGetOrder(SNESLineSearch linesearch,PetscInt * order)1293d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch linesearch, PetscInt *order)
1294d71ae5a4SJacob Faibussowitsch {
129559405d5eSPeter Brune PetscFunctionBegin;
129659405d5eSPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
12974f572ea9SToby Isaac PetscAssertPointer(order, 2);
129859405d5eSPeter Brune *order = linesearch->order;
12993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
130059405d5eSPeter Brune }
130159405d5eSPeter Brune
130259405d5eSPeter Brune /*@
13031f8196a2SJed Brown SNESLineSearchSetOrder - Sets the maximum order of the polynomial fit used in the line search
130459405d5eSPeter Brune
130559405d5eSPeter Brune Input Parameters:
1306420bcc1bSBarry Smith + linesearch - the line search context
1307ceaaa498SBarry Smith - order - The order
130859405d5eSPeter Brune
130959405d5eSPeter Brune Level: intermediate
131059405d5eSPeter Brune
1311ceaaa498SBarry Smith Values for `order`\:
1312f6dfbefdSBarry Smith + 1 or `SNES_LINESEARCH_ORDER_LINEAR` - linear order
1313f6dfbefdSBarry Smith . 2 or `SNES_LINESEARCH_ORDER_QUADRATIC` - quadratic order
1314f6dfbefdSBarry Smith - 3 or `SNES_LINESEARCH_ORDER_CUBIC` - cubic order
131578bcb3b5SPeter Brune
1316420bcc1bSBarry Smith Options Database Key:
1317420bcc1bSBarry Smith . -snes_linesearch_order <order> - 1, 2, 3. Most types only support certain orders (`SNESLINESEARCHBT` supports 2 or 3)
1318420bcc1bSBarry Smith
1319ceaaa498SBarry Smith Note:
1320ceaaa498SBarry Smith These orders are supported by `SNESLINESEARCHBT` and `SNESLINESEARCHCP`
132159405d5eSPeter Brune
1322420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetOrder()`, `SNESLineSearchSetDamping()`
132359405d5eSPeter Brune @*/
SNESLineSearchSetOrder(SNESLineSearch linesearch,PetscInt order)1324d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch linesearch, PetscInt order)
1325d71ae5a4SJacob Faibussowitsch {
132659405d5eSPeter Brune PetscFunctionBegin;
132759405d5eSPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
132859405d5eSPeter Brune linesearch->order = order;
13293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
133059405d5eSPeter Brune }
133159405d5eSPeter Brune
1332f40b411bSPeter Brune /*@
1333420bcc1bSBarry Smith SNESLineSearchGetNorms - Gets the norms for the current solution `X`, the current update `Y`, and the current function value `F`.
1334f40b411bSPeter Brune
1335f6dfbefdSBarry Smith Not Collective
1336f6dfbefdSBarry Smith
1337f899ff85SJose E. Roman Input Parameter:
1338420bcc1bSBarry Smith . linesearch - the line search context
1339f40b411bSPeter Brune
1340f40b411bSPeter Brune Output Parameters:
1341f40b411bSPeter Brune + xnorm - The norm of the current solution
1342a68bbae5SBarry Smith . fnorm - The norm of the current function, this is the `norm(function(X))` where `X` is the current solution.
1343a99ef635SJonas Heinzmann - ynorm - The norm of the current update (after scaling by the linesearch computed `lambda`)
1344f40b411bSPeter Brune
134578bcb3b5SPeter Brune Level: developer
1346f40b411bSPeter Brune
1347420bcc1bSBarry Smith Notes:
1348420bcc1bSBarry Smith Some values may not be up-to-date at particular points in the code.
1349a68bbae5SBarry Smith
1350a68bbae5SBarry Smith This, in combination with `SNESLineSearchSetNorms()`, allow the line search and the `SNESSolve_XXX()` to share
1351a68bbae5SBarry Smith computed values.
1352a68bbae5SBarry Smith
13530241b274SPierre Jolivet .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetNorms()`, `SNESLineSearchGetVecs()`
1354f40b411bSPeter Brune @*/
SNESLineSearchGetNorms(SNESLineSearch linesearch,PetscReal * xnorm,PetscReal * fnorm,PetscReal * ynorm)1355d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal *xnorm, PetscReal *fnorm, PetscReal *ynorm)
1356d71ae5a4SJacob Faibussowitsch {
1357bf7f4e0aSPeter Brune PetscFunctionBegin;
1358f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
1359f5af7f23SKarl Rupp if (xnorm) *xnorm = linesearch->xnorm;
1360f5af7f23SKarl Rupp if (fnorm) *fnorm = linesearch->fnorm;
1361f5af7f23SKarl Rupp if (ynorm) *ynorm = linesearch->ynorm;
13623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1363bf7f4e0aSPeter Brune }
1364bf7f4e0aSPeter Brune
1365f40b411bSPeter Brune /*@
13661bd63e3eSSatish Balay SNESLineSearchSetNorms - Sets the computed norms for the current solution `X`, the current update `Y`, and the current function value `F`.
1367f40b411bSPeter Brune
1368c3339decSBarry Smith Collective
1369f6dfbefdSBarry Smith
1370f40b411bSPeter Brune Input Parameters:
1371420bcc1bSBarry Smith + linesearch - the line search context
1372f40b411bSPeter Brune . xnorm - The norm of the current solution
1373a68bbae5SBarry Smith . fnorm - The norm of the current function, this is the `norm(function(X))` where `X` is the current solution
1374a99ef635SJonas Heinzmann - ynorm - The norm of the current update (after scaling by the linesearch computed `lambda`)
1375f40b411bSPeter Brune
1376f6dfbefdSBarry Smith Level: developer
1377f40b411bSPeter Brune
1378420bcc1bSBarry Smith Note:
1379420bcc1bSBarry Smith This is called by the line search routines to store the values they have just computed
1380420bcc1bSBarry Smith
1381420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()`
1382f40b411bSPeter Brune @*/
SNESLineSearchSetNorms(SNESLineSearch linesearch,PetscReal xnorm,PetscReal fnorm,PetscReal ynorm)1383d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
1384d71ae5a4SJacob Faibussowitsch {
13856a388c36SPeter Brune PetscFunctionBegin;
1386f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
13876a388c36SPeter Brune linesearch->xnorm = xnorm;
13886a388c36SPeter Brune linesearch->fnorm = fnorm;
13896a388c36SPeter Brune linesearch->ynorm = ynorm;
13903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
13916a388c36SPeter Brune }
13926a388c36SPeter Brune
1393f40b411bSPeter Brune /*@
1394420bcc1bSBarry Smith SNESLineSearchComputeNorms - Explicitly computes the norms of the current solution `X`, the current update `Y`, and the current function value `F`.
1395f40b411bSPeter Brune
13962fe279fdSBarry Smith Input Parameter:
1397420bcc1bSBarry Smith . linesearch - the line search context
1398f40b411bSPeter Brune
139920f4b53cSBarry Smith Options Database Key:
14008e557f58SPeter Brune . -snes_linesearch_norms - turn norm computation on or off
1401f40b411bSPeter Brune
1402f40b411bSPeter Brune Level: intermediate
1403f40b411bSPeter Brune
1404420bcc1bSBarry Smith Developer Note:
1405420bcc1bSBarry Smith The options database key is misnamed. It should be -snes_linesearch_compute_norms
1406420bcc1bSBarry Smith
1407420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms`, `SNESLineSearchSetNorms()`, `SNESLineSearchSetComputeNorms()`
1408f40b411bSPeter Brune @*/
SNESLineSearchComputeNorms(SNESLineSearch linesearch)1409d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
1410d71ae5a4SJacob Faibussowitsch {
14119bd66eb0SPeter Brune SNES snes;
1412bf388a1fSBarry Smith
14136a388c36SPeter Brune PetscFunctionBegin;
14146a388c36SPeter Brune if (linesearch->norms) {
14159bd66eb0SPeter Brune if (linesearch->ops->vinorm) {
14169566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
14179566063dSJacob Faibussowitsch PetscCall(VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
14189566063dSJacob Faibussowitsch PetscCall(VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm));
14199566063dSJacob Faibussowitsch PetscCall((*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm));
14209bd66eb0SPeter Brune } else {
14219566063dSJacob Faibussowitsch PetscCall(VecNormBegin(linesearch->vec_func, NORM_2, &linesearch->fnorm));
14229566063dSJacob Faibussowitsch PetscCall(VecNormBegin(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
14239566063dSJacob Faibussowitsch PetscCall(VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm));
14249566063dSJacob Faibussowitsch PetscCall(VecNormEnd(linesearch->vec_func, NORM_2, &linesearch->fnorm));
14259566063dSJacob Faibussowitsch PetscCall(VecNormEnd(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
14269566063dSJacob Faibussowitsch PetscCall(VecNormEnd(linesearch->vec_update, NORM_2, &linesearch->ynorm));
14276a388c36SPeter Brune }
14289bd66eb0SPeter Brune }
14293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
14306a388c36SPeter Brune }
14316a388c36SPeter Brune
14326f263ca3SPeter Brune /*@
14336f263ca3SPeter Brune SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search.
14346f263ca3SPeter Brune
14356f263ca3SPeter Brune Input Parameters:
1436420bcc1bSBarry Smith + linesearch - the line search context
143778bcb3b5SPeter Brune - flg - indicates whether or not to compute norms
14386f263ca3SPeter Brune
143920f4b53cSBarry Smith Options Database Key:
1440420bcc1bSBarry Smith . -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic (none) `SNESLINESEARCHBASIC` line search
14416f263ca3SPeter Brune
144220f4b53cSBarry Smith Level: intermediate
144320f4b53cSBarry Smith
1444f6dfbefdSBarry Smith Note:
1445f6dfbefdSBarry 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.
14466f263ca3SPeter Brune
1447420bcc1bSBarry Smith Developer Note:
1448420bcc1bSBarry Smith The options database key is misnamed. It should be -snes_linesearch_compute_norms
1449420bcc1bSBarry Smith
1450420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetNorms()`, `SNESLineSearchComputeNorms()`, `SNESLINESEARCHBASIC`
14516f263ca3SPeter Brune @*/
SNESLineSearchSetComputeNorms(SNESLineSearch linesearch,PetscBool flg)1452d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg)
1453d71ae5a4SJacob Faibussowitsch {
14546f263ca3SPeter Brune PetscFunctionBegin;
14556f263ca3SPeter Brune linesearch->norms = flg;
14563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
14576f263ca3SPeter Brune }
14586f263ca3SPeter Brune
1459f40b411bSPeter Brune /*@
1460f6dfbefdSBarry Smith SNESLineSearchGetVecs - Gets the vectors from the `SNESLineSearch` context
1461f6dfbefdSBarry Smith
14628f14a041SBarry Smith Not Collective but the vectors are parallel
1463f40b411bSPeter Brune
1464f899ff85SJose E. Roman Input Parameter:
1465420bcc1bSBarry Smith . linesearch - the line search context
1466f40b411bSPeter Brune
1467f40b411bSPeter Brune Output Parameters:
14686232e825SPeter Brune + X - Solution vector
14696232e825SPeter Brune . F - Function vector
14706232e825SPeter Brune . Y - Search direction vector
14716232e825SPeter Brune . W - Solution work vector
14726232e825SPeter Brune - G - Function work vector
14736232e825SPeter Brune
147420f4b53cSBarry Smith Level: advanced
147520f4b53cSBarry Smith
14767bba9028SPeter Brune Notes:
147720f4b53cSBarry Smith At the beginning of a line search application, `X` should contain a
147820f4b53cSBarry Smith solution and the vector `F` the function computed at `X`. At the end of the
147920f4b53cSBarry Smith line search application, `X` should contain the new solution, and `F` the
14806232e825SPeter Brune function evaluated at the new solution.
1481f40b411bSPeter Brune
1482f6dfbefdSBarry Smith These vectors are owned by the `SNESLineSearch` and should not be destroyed by the caller
14832a7a6963SBarry Smith
1484420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()`
1485f40b411bSPeter Brune @*/
SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec * X,Vec * F,Vec * Y,Vec * W,Vec * G)1486d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch, Vec *X, Vec *F, Vec *Y, Vec *W, Vec *G)
1487d71ae5a4SJacob Faibussowitsch {
14886a388c36SPeter Brune PetscFunctionBegin;
1489f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
14906a388c36SPeter Brune if (X) {
14914f572ea9SToby Isaac PetscAssertPointer(X, 2);
14926a388c36SPeter Brune *X = linesearch->vec_sol;
14936a388c36SPeter Brune }
14946a388c36SPeter Brune if (F) {
14954f572ea9SToby Isaac PetscAssertPointer(F, 3);
14966a388c36SPeter Brune *F = linesearch->vec_func;
14976a388c36SPeter Brune }
14986a388c36SPeter Brune if (Y) {
14994f572ea9SToby Isaac PetscAssertPointer(Y, 4);
15006a388c36SPeter Brune *Y = linesearch->vec_update;
15016a388c36SPeter Brune }
15026a388c36SPeter Brune if (W) {
15034f572ea9SToby Isaac PetscAssertPointer(W, 5);
15046a388c36SPeter Brune *W = linesearch->vec_sol_new;
15056a388c36SPeter Brune }
15066a388c36SPeter Brune if (G) {
15074f572ea9SToby Isaac PetscAssertPointer(G, 6);
15086a388c36SPeter Brune *G = linesearch->vec_func_new;
15096a388c36SPeter Brune }
15103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
15116a388c36SPeter Brune }
15126a388c36SPeter Brune
1513f40b411bSPeter Brune /*@
1514f6dfbefdSBarry Smith SNESLineSearchSetVecs - Sets the vectors on the `SNESLineSearch` context
1515f6dfbefdSBarry Smith
1516c3339decSBarry Smith Logically Collective
1517f40b411bSPeter Brune
1518f40b411bSPeter Brune Input Parameters:
1519420bcc1bSBarry Smith + linesearch - the line search context
15206232e825SPeter Brune . X - Solution vector
15216232e825SPeter Brune . F - Function vector
15226232e825SPeter Brune . Y - Search direction vector
15236232e825SPeter Brune . W - Solution work vector
15246232e825SPeter Brune - G - Function work vector
1525f40b411bSPeter Brune
1526420bcc1bSBarry Smith Level: developer
1527f40b411bSPeter Brune
1528420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetNorms()`, `SNESLineSearchGetVecs()`
1529f40b411bSPeter Brune @*/
SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W,Vec G)1530d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch, Vec X, Vec F, Vec Y, Vec W, Vec G)
1531d71ae5a4SJacob Faibussowitsch {
15326a388c36SPeter Brune PetscFunctionBegin;
1533f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
15346a388c36SPeter Brune if (X) {
15356a388c36SPeter Brune PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
15366a388c36SPeter Brune linesearch->vec_sol = X;
15376a388c36SPeter Brune }
15386a388c36SPeter Brune if (F) {
15396a388c36SPeter Brune PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
15406a388c36SPeter Brune linesearch->vec_func = F;
15416a388c36SPeter Brune }
15426a388c36SPeter Brune if (Y) {
15436a388c36SPeter Brune PetscValidHeaderSpecific(Y, VEC_CLASSID, 4);
15446a388c36SPeter Brune linesearch->vec_update = Y;
15456a388c36SPeter Brune }
15466a388c36SPeter Brune if (W) {
15476a388c36SPeter Brune PetscValidHeaderSpecific(W, VEC_CLASSID, 5);
15486a388c36SPeter Brune linesearch->vec_sol_new = W;
15496a388c36SPeter Brune }
15506a388c36SPeter Brune if (G) {
15516a388c36SPeter Brune PetscValidHeaderSpecific(G, VEC_CLASSID, 6);
15526a388c36SPeter Brune linesearch->vec_func_new = G;
15536a388c36SPeter Brune }
15543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
15556a388c36SPeter Brune }
15566a388c36SPeter Brune
1557cc4c1da9SBarry Smith /*@
1558f1c6b773SPeter Brune SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1559f6dfbefdSBarry Smith `SNESLineSearch` options in the database.
1560e7058c64SPeter Brune
1561c3339decSBarry Smith Logically Collective
1562e7058c64SPeter Brune
1563e7058c64SPeter Brune Input Parameters:
1564f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context
1565e7058c64SPeter Brune - prefix - the prefix to prepend to all option names
1566e7058c64SPeter Brune
156720f4b53cSBarry Smith Level: advanced
156820f4b53cSBarry Smith
1569f6dfbefdSBarry Smith Note:
1570e7058c64SPeter Brune A hyphen (-) must NOT be given at the beginning of the prefix name.
1571e7058c64SPeter Brune The first character of all runtime options is AUTOMATICALLY the hyphen.
1572e7058c64SPeter Brune
1573420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch()`, `SNESLineSearchSetFromOptions()`, `SNESGetOptionsPrefix()`
1574e7058c64SPeter Brune @*/
SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[])1575d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch, const char prefix[])
1576d71ae5a4SJacob Faibussowitsch {
1577e7058c64SPeter Brune PetscFunctionBegin;
1578f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
15799566063dSJacob Faibussowitsch PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)linesearch, prefix));
15803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1581e7058c64SPeter Brune }
1582e7058c64SPeter Brune
1583cc4c1da9SBarry Smith /*@
1584f6dfbefdSBarry Smith SNESLineSearchGetOptionsPrefix - Gets the prefix used for searching for all
1585f1c6b773SPeter Brune SNESLineSearch options in the database.
1586e7058c64SPeter Brune
1587e7058c64SPeter Brune Not Collective
1588e7058c64SPeter Brune
1589e7058c64SPeter Brune Input Parameter:
1590f6dfbefdSBarry Smith . linesearch - the `SNESLineSearch` context
1591e7058c64SPeter Brune
1592e7058c64SPeter Brune Output Parameter:
1593e7058c64SPeter Brune . prefix - pointer to the prefix string used
1594e7058c64SPeter Brune
1595e7058c64SPeter Brune Level: advanced
1596e7058c64SPeter Brune
1597420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESAppendOptionsPrefix()`
1598e7058c64SPeter Brune @*/
SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char * prefix[])1599d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch, const char *prefix[])
1600d71ae5a4SJacob Faibussowitsch {
1601e7058c64SPeter Brune PetscFunctionBegin;
1602f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
16039566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)linesearch, prefix));
16043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1605e7058c64SPeter Brune }
1606bf7f4e0aSPeter Brune
16078d359177SBarry Smith /*@C
1608f6dfbefdSBarry Smith SNESLineSearchSetWorkVecs - Sets work vectors for the line search.
1609f40b411bSPeter Brune
1610d8d19677SJose E. Roman Input Parameters:
1611f6dfbefdSBarry Smith + linesearch - the `SNESLineSearch` context
1612f40b411bSPeter Brune - nwork - the number of work vectors
1613f40b411bSPeter Brune
1614f40b411bSPeter Brune Level: developer
1615f40b411bSPeter Brune
1616420bcc1bSBarry Smith Developer Note:
1617420bcc1bSBarry Smith This is called from within the set up routines for each of the line search types `SNESLineSearchType`
1618420bcc1bSBarry Smith
1619420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESSetWorkVecs()`
1620f40b411bSPeter Brune @*/
SNESLineSearchSetWorkVecs(SNESLineSearch linesearch,PetscInt nwork)1621d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork)
1622d71ae5a4SJacob Faibussowitsch {
1623bf7f4e0aSPeter Brune PetscFunctionBegin;
16240fdf79fbSJacob Faibussowitsch PetscCheck(linesearch->vec_sol, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
16259566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work));
16263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1627bf7f4e0aSPeter Brune }
1628bf7f4e0aSPeter Brune
1629f40b411bSPeter Brune /*@
1630422a814eSBarry Smith SNESLineSearchGetReason - Gets the success/failure status of the last line search application
1631f40b411bSPeter Brune
16322fe279fdSBarry Smith Input Parameter:
1633420bcc1bSBarry Smith . linesearch - the line search context
1634f40b411bSPeter Brune
16352fe279fdSBarry Smith Output Parameter:
163676c63389SBarry Smith . reason - The success or failure status
1637f40b411bSPeter Brune
163820f4b53cSBarry Smith Level: developer
163920f4b53cSBarry Smith
1640f6dfbefdSBarry Smith Note:
1641420bcc1bSBarry Smith This is typically called after `SNESLineSearchApply()` in order to determine if the line search failed
1642420bcc1bSBarry Smith (and set into the `SNES` convergence accordingly).
1643cd7522eaSPeter Brune
1644420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetReason()`, `SNESLineSearchReason`
1645f40b411bSPeter Brune @*/
SNESLineSearchGetReason(SNESLineSearch linesearch,SNESLineSearchReason * reason)164676c63389SBarry Smith PetscErrorCode SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *reason)
1647d71ae5a4SJacob Faibussowitsch {
1648bf7f4e0aSPeter Brune PetscFunctionBegin;
1649f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
165076c63389SBarry Smith PetscAssertPointer(reason, 2);
165176c63389SBarry Smith *reason = linesearch->reason;
16523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1653bf7f4e0aSPeter Brune }
1654bf7f4e0aSPeter Brune
1655f40b411bSPeter Brune /*@
165676c63389SBarry Smith SNESLineSearchSetReason - Sets the success/failure reason of the line search application
1657f40b411bSPeter Brune
16585d83a8b1SBarry Smith Logically Collective; No Fortran Support
16595d83a8b1SBarry Smith
1660f40b411bSPeter Brune Input Parameters:
1661420bcc1bSBarry Smith + linesearch - the line search context
166276c63389SBarry Smith - reason - The success or failure reason
1663f40b411bSPeter Brune
166420f4b53cSBarry Smith Level: developer
166520f4b53cSBarry Smith
166676c63389SBarry Smith Notes:
1667420bcc1bSBarry Smith This is typically called in a `SNESLineSearchType` implementation of `SNESLineSearchApply()` or a `SNESLINESEARCHSHELL` implementation to set
1668cd7522eaSPeter Brune the success or failure of the line search method.
1669cd7522eaSPeter Brune
167076c63389SBarry Smith Do not call this from callbacks provided with `SNESSetFunction()`, instead perhaps use `SNESSetFunctionDomainError()`
167176c63389SBarry Smith
167276c63389SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchReason`, `SNESLineSearchGetSResult()`, `SNESSetFunctionDomainError()`, `SNESSetFunction()`
1673f40b411bSPeter Brune @*/
SNESLineSearchSetReason(SNESLineSearch linesearch,SNESLineSearchReason reason)167476c63389SBarry Smith PetscErrorCode SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason reason)
1675d71ae5a4SJacob Faibussowitsch {
16766a388c36SPeter Brune PetscFunctionBegin;
16775fd66863SKarl Rupp PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
167876c63389SBarry Smith linesearch->reason = reason;
16793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
16806a388c36SPeter Brune }
16816a388c36SPeter Brune
1682ceaaa498SBarry Smith // PetscClangLinter pragma disable: -fdoc-param-list-func-parameter-documentation
16839bd66eb0SPeter Brune /*@C
1684f1c6b773SPeter Brune SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.
16859bd66eb0SPeter Brune
1686c3339decSBarry Smith Logically Collective
1687f6dfbefdSBarry Smith
16889bd66eb0SPeter Brune Input Parameters:
1689ceaaa498SBarry Smith + linesearch - the linesearch object
16908434afd1SBarry Smith . projectfunc - function for projecting the function to the bounds, see `SNESLineSearchVIProjectFn` for calling sequence
1691d5def619SJonas Heinzmann . normfunc - function for computing the norm of an active set, see `SNESLineSearchVINormFn` for calling sequence
1692d5def619SJonas Heinzmann - dirderivfunc - function for computing the directional derivative of an active set, see `SNESLineSearchVIDirDerivFn` for calling sequence
16939bd66eb0SPeter Brune
1694f6dfbefdSBarry Smith Level: advanced
16959bd66eb0SPeter Brune
1696ceaaa498SBarry Smith Notes:
169720f4b53cSBarry Smith The VI solvers require projection of the solution to the feasible set. `projectfunc` should implement this.
169820f4b53cSBarry Smith
169920f4b53cSBarry Smith The VI solvers require special evaluation of the function norm such that the norm is only calculated
170020f4b53cSBarry Smith on the inactive set. This should be implemented by `normfunc`.
170120f4b53cSBarry Smith
1702d5def619SJonas Heinzmann The VI solvers further require special evaluation of the directional derivative (when assuming that there exists some $G(x)$
1703d5def619SJonas Heinzmann for which the `SNESFunctionFn` $F(x) = grad G(x)$) such that it is only calculated on the inactive set.
1704d5def619SJonas Heinzmann This should be implemented by `dirderivfunc`.
1705d5def619SJonas Heinzmann
17069bcc50f1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetVIFunctions()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()`,
1707d5def619SJonas Heinzmann `SNESLineSearchVIProjectFn`, `SNESLineSearchVINormFn`, `SNESLineSearchVIDirDerivFn`
17089bd66eb0SPeter Brune @*/
SNESLineSearchSetVIFunctions(SNESLineSearch linesearch,SNESLineSearchVIProjectFn * projectfunc,SNESLineSearchVINormFn * normfunc,SNESLineSearchVIDirDerivFn * dirderivfunc)1709d5def619SJonas Heinzmann PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFn *projectfunc, SNESLineSearchVINormFn *normfunc, SNESLineSearchVIDirDerivFn *dirderivfunc)
1710d71ae5a4SJacob Faibussowitsch {
17119bd66eb0SPeter Brune PetscFunctionBegin;
1712f1c6b773SPeter Brune PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
17139bd66eb0SPeter Brune if (projectfunc) linesearch->ops->viproject = projectfunc;
17149bd66eb0SPeter Brune if (normfunc) linesearch->ops->vinorm = normfunc;
1715d5def619SJonas Heinzmann if (dirderivfunc) linesearch->ops->vidirderiv = dirderivfunc;
17163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
17179bd66eb0SPeter Brune }
17189bd66eb0SPeter Brune
17199bd66eb0SPeter Brune /*@C
1720f1c6b773SPeter Brune SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.
17219bd66eb0SPeter Brune
1722f6dfbefdSBarry Smith Not Collective
1723f6dfbefdSBarry Smith
1724f899ff85SJose E. Roman Input Parameter:
1725f6dfbefdSBarry Smith . linesearch - the line search context, obtain with `SNESGetLineSearch()`
17269bd66eb0SPeter Brune
17279bd66eb0SPeter Brune Output Parameters:
17288434afd1SBarry Smith + projectfunc - function for projecting the function to the bounds, see `SNESLineSearchVIProjectFn` for calling sequence
1729d5def619SJonas Heinzmann . normfunc - function for computing the norm of an active set, see `SNESLineSearchVINormFn ` for calling sequence
1730d5def619SJonas Heinzmann - dirderivfunc - function for computing the directional derivative of an active set, see `SNESLineSearchVIDirDerivFn` for calling sequence
17319bd66eb0SPeter Brune
1732f6dfbefdSBarry Smith Level: advanced
17339bd66eb0SPeter Brune
17349bcc50f1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetVIFunctions()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()`,
17358434afd1SBarry Smith `SNESLineSearchVIProjectFn`, `SNESLineSearchVINormFn`
17369bd66eb0SPeter Brune @*/
SNESLineSearchGetVIFunctions(SNESLineSearch linesearch,SNESLineSearchVIProjectFn ** projectfunc,SNESLineSearchVINormFn ** normfunc,SNESLineSearchVIDirDerivFn ** dirderivfunc)1737d5def619SJonas Heinzmann PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFn **projectfunc, SNESLineSearchVINormFn **normfunc, SNESLineSearchVIDirDerivFn **dirderivfunc)
1738d71ae5a4SJacob Faibussowitsch {
17399bd66eb0SPeter Brune PetscFunctionBegin;
17409bd66eb0SPeter Brune if (projectfunc) *projectfunc = linesearch->ops->viproject;
17419bd66eb0SPeter Brune if (normfunc) *normfunc = linesearch->ops->vinorm;
1742d5def619SJonas Heinzmann if (dirderivfunc) *dirderivfunc = linesearch->ops->vidirderiv;
17433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
17449bd66eb0SPeter Brune }
17459bd66eb0SPeter Brune
1746bf7f4e0aSPeter Brune /*@C
1747420bcc1bSBarry Smith SNESLineSearchRegister - register a line search type `SNESLineSearchType`
1748ceaaa498SBarry Smith
1749cc4c1da9SBarry Smith Logically Collective, No Fortran Support
1750cc4c1da9SBarry Smith
1751ceaaa498SBarry Smith Input Parameters:
1752420bcc1bSBarry Smith + sname - name of the `SNESLineSearchType()`
1753ceaaa498SBarry Smith - function - the creation function for that type
1754ceaaa498SBarry Smith
1755ceaaa498SBarry Smith Calling sequence of `function`:
1756ceaaa498SBarry Smith . ls - the line search context
1757bf7f4e0aSPeter Brune
1758bf7f4e0aSPeter Brune Level: advanced
1759f6dfbefdSBarry Smith
1760420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchSetType()`
1761bf7f4e0aSPeter Brune @*/
SNESLineSearchRegister(const char sname[],PetscErrorCode (* function)(SNESLineSearch ls))1762ceaaa498SBarry Smith PetscErrorCode SNESLineSearchRegister(const char sname[], PetscErrorCode (*function)(SNESLineSearch ls))
1763d71ae5a4SJacob Faibussowitsch {
1764bf7f4e0aSPeter Brune PetscFunctionBegin;
17659566063dSJacob Faibussowitsch PetscCall(SNESInitializePackage());
17669566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&SNESLineSearchList, sname, function));
17673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1768bf7f4e0aSPeter Brune }
1769