Actual source code: linesearch.c
1: #include <petsc/private/linesearchimpl.h>
3: PetscBool SNESLineSearchRegisterAllCalled = PETSC_FALSE;
4: PetscFunctionList SNESLineSearchList = NULL;
6: PetscClassId SNESLINESEARCH_CLASSID;
7: PetscLogEvent SNESLINESEARCH_Apply;
9: /*@
10: SNESLineSearchMonitorCancel - Clears all the monitor functions for a `SNESLineSearch` object.
12: Logically Collective
14: Input Parameter:
15: . ls - the `SNESLineSearch` context
17: Options Database Key:
18: . -snes_linesearch_monitor_cancel - cancels all monitors that have been hardwired
19: into a code by calls to `SNESLineSearchMonitorSet()`, but does not cancel those
20: set via the options database
22: Level: advanced
24: Notes:
25: There is no way to clear one specific monitor from a `SNESLineSearch` object.
27: This does not clear the monitor set with `SNESLineSearchSetDefaultMonitor()` use `SNESLineSearchSetDefaultMonitor`(ls,NULL) to cancel
28: that one.
30: .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorDefault()`, `SNESLineSearchMonitorSet()`
31: @*/
32: PetscErrorCode SNESLineSearchMonitorCancel(SNESLineSearch ls)
33: {
34: PetscInt i;
36: PetscFunctionBegin;
38: for (i = 0; i < ls->numbermonitors; i++) {
39: if (ls->monitordestroy[i]) PetscCall((*ls->monitordestroy[i])(&ls->monitorcontext[i]));
40: }
41: ls->numbermonitors = 0;
42: PetscFunctionReturn(PETSC_SUCCESS);
43: }
45: /*@
46: SNESLineSearchMonitor - runs the user provided monitor routines, if they exist
48: Collective
50: Input Parameter:
51: . ls - the linesearch object
53: Level: developer
55: Note:
56: This routine is called by the `SNES` implementations.
57: It does not typically need to be called by the user.
59: .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorSet()`
60: @*/
61: PetscErrorCode SNESLineSearchMonitor(SNESLineSearch ls)
62: {
63: PetscInt i, n = ls->numbermonitors;
65: PetscFunctionBegin;
66: for (i = 0; i < n; i++) PetscCall((*ls->monitorftns[i])(ls, ls->monitorcontext[i]));
67: PetscFunctionReturn(PETSC_SUCCESS);
68: }
70: /*@C
71: SNESLineSearchMonitorSet - Sets an ADDITIONAL function that is to be used at every
72: iteration of the nonlinear solver to display the iteration's
73: progress.
75: Logically Collective
77: Input Parameters:
78: + ls - the `SNESLineSearch` context
79: . f - the monitor function
80: . mctx - [optional] user-defined context for private data for the
81: monitor routine (use `NULL` if no context is desired)
82: - monitordestroy - [optional] routine that frees monitor context
83: (may be `NULL`)
85: Level: intermediate
87: Note:
88: Several different monitoring routines may be set by calling
89: `SNESLineSearchMonitorSet()` multiple times; all will be called in the
90: order in which they were set.
92: Fortran Note:
93: Only a single monitor function can be set for each `SNESLineSearch` object
95: .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorDefault()`, `SNESLineSearchMonitorCancel()`
96: @*/
97: PetscErrorCode SNESLineSearchMonitorSet(SNESLineSearch ls, PetscErrorCode (*f)(SNESLineSearch, void *), void *mctx, PetscErrorCode (*monitordestroy)(void **))
98: {
99: PetscInt i;
100: PetscBool identical;
102: PetscFunctionBegin;
104: for (i = 0; i < ls->numbermonitors; i++) {
105: PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))f, mctx, monitordestroy, (PetscErrorCode(*)(void))ls->monitorftns[i], ls->monitorcontext[i], ls->monitordestroy[i], &identical));
106: if (identical) PetscFunctionReturn(PETSC_SUCCESS);
107: }
108: PetscCheck(ls->numbermonitors < MAXSNESLSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many monitors set");
109: ls->monitorftns[ls->numbermonitors] = f;
110: ls->monitordestroy[ls->numbermonitors] = monitordestroy;
111: ls->monitorcontext[ls->numbermonitors++] = (void *)mctx;
112: PetscFunctionReturn(PETSC_SUCCESS);
113: }
115: /*@C
116: SNESLineSearchMonitorSolutionUpdate - Monitors each update of the function value the linesearch tries
118: Collective
120: Input Parameters:
121: + ls - the `SNES` linesearch object
122: - vf - the context for the monitor, in this case it is an `PetscViewerAndFormat`
124: Level: intermediate
126: .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESMonitorSet()`, `SNESMonitorSolution()`
127: @*/
128: PetscErrorCode SNESLineSearchMonitorSolutionUpdate(SNESLineSearch ls, PetscViewerAndFormat *vf)
129: {
130: PetscViewer viewer = vf->viewer;
131: Vec Y, W, G;
133: PetscFunctionBegin;
134: PetscCall(SNESLineSearchGetVecs(ls, NULL, NULL, &Y, &W, &G));
135: PetscCall(PetscViewerPushFormat(viewer, vf->format));
136: PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted update to solution \n"));
137: PetscCall(VecView(Y, viewer));
138: PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted new solution \n"));
139: PetscCall(VecView(W, viewer));
140: PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted updated function value\n"));
141: PetscCall(VecView(G, viewer));
142: PetscCall(PetscViewerPopFormat(viewer));
143: PetscFunctionReturn(PETSC_SUCCESS);
144: }
146: /*@
147: SNESLineSearchCreate - Creates the line search context.
149: Logically Collective
151: Input Parameter:
152: . comm - MPI communicator for the line search (typically from the associated `SNES` context).
154: Output Parameter:
155: . outlinesearch - the new linesearch context
157: Level: developer
159: Note:
160: The preferred calling sequence for users is to use `SNESGetLineSearch()` to acquire the `SNESLineSearch` instance
161: already associated with the `SNES`.
163: .seealso: `SNESLineSearch`, `LineSearchDestroy()`, `SNESGetLineSearch()`
164: @*/
165: PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch)
166: {
167: SNESLineSearch linesearch;
169: PetscFunctionBegin;
171: PetscCall(SNESInitializePackage());
172: *outlinesearch = NULL;
174: PetscCall(PetscHeaderCreate(linesearch, SNESLINESEARCH_CLASSID, "SNESLineSearch", "Linesearch", "SNESLineSearch", comm, SNESLineSearchDestroy, SNESLineSearchView));
176: linesearch->vec_sol_new = NULL;
177: linesearch->vec_func_new = NULL;
178: linesearch->vec_sol = NULL;
179: linesearch->vec_func = NULL;
180: linesearch->vec_update = NULL;
182: linesearch->lambda = 1.0;
183: linesearch->fnorm = 1.0;
184: linesearch->ynorm = 1.0;
185: linesearch->xnorm = 1.0;
186: linesearch->result = SNES_LINESEARCH_SUCCEEDED;
187: linesearch->norms = PETSC_TRUE;
188: linesearch->keeplambda = PETSC_FALSE;
189: linesearch->damping = 1.0;
190: linesearch->maxstep = 1e8;
191: linesearch->steptol = 1e-12;
192: linesearch->rtol = 1e-8;
193: linesearch->atol = 1e-15;
194: linesearch->ltol = 1e-8;
195: linesearch->precheckctx = NULL;
196: linesearch->postcheckctx = NULL;
197: linesearch->max_its = 1;
198: linesearch->setupcalled = PETSC_FALSE;
199: linesearch->monitor = NULL;
200: *outlinesearch = linesearch;
201: PetscFunctionReturn(PETSC_SUCCESS);
202: }
204: /*@
205: SNESLineSearchSetUp - Prepares the line search for being applied by allocating
206: any required vectors.
208: Collective
210: Input Parameter:
211: . linesearch - The `SNESLineSearch` instance.
213: Level: advanced
215: Note:
216: For most cases, this needn't be called by users or outside of `SNESLineSearchApply()`.
217: The only current case where this is called outside of this is for the VI
218: solvers, which modify the solution and work vectors before the first call
219: of `SNESLineSearchApply()`, requiring the `SNESLineSearch` work vectors to be
220: allocated upfront.
222: .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchReset()`
223: @*/
225: PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch)
226: {
227: PetscFunctionBegin;
228: if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC));
229: if (!linesearch->setupcalled) {
230: if (!linesearch->vec_sol_new) PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new));
231: if (!linesearch->vec_func_new) PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new));
232: if (linesearch->ops->setup) PetscUseTypeMethod(linesearch, setup);
233: if (!linesearch->ops->snesfunc) PetscCall(SNESLineSearchSetFunction(linesearch, SNESComputeFunction));
234: linesearch->lambda = linesearch->damping;
235: linesearch->setupcalled = PETSC_TRUE;
236: }
237: PetscFunctionReturn(PETSC_SUCCESS);
238: }
240: /*@
241: SNESLineSearchReset - Undoes the `SNESLineSearchSetUp()` and deletes any `Vec`s or `Mat`s allocated by the line search.
243: Collective
245: Input Parameter:
246: . linesearch - The `SNESLineSearch` instance.
248: Level: developer
250: Note:
251: Usually only called by `SNESReset()`
253: .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchSetUp()`
254: @*/
256: PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch)
257: {
258: PetscFunctionBegin;
259: if (linesearch->ops->reset) PetscUseTypeMethod(linesearch, reset);
261: PetscCall(VecDestroy(&linesearch->vec_sol_new));
262: PetscCall(VecDestroy(&linesearch->vec_func_new));
264: PetscCall(VecDestroyVecs(linesearch->nwork, &linesearch->work));
266: linesearch->nwork = 0;
267: linesearch->setupcalled = PETSC_FALSE;
268: PetscFunctionReturn(PETSC_SUCCESS);
269: }
271: /*@C
272: SNESLineSearchSetFunction - Sets the function evaluation used by the `SNES` line search
273: `
274: Input Parameters:
275: . linesearch - the `SNESLineSearch` context
276: + func - function evaluation routine, this is usually the function provided with `SNESSetFunction()`
278: Level: developer
280: .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESSetFunction()`
281: @*/
282: PetscErrorCode SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES, Vec, Vec))
283: {
284: PetscFunctionBegin;
286: linesearch->ops->snesfunc = func;
287: PetscFunctionReturn(PETSC_SUCCESS);
288: }
290: /*@C
291: SNESLineSearchSetPreCheck - Sets a user function that is called after the initial search direction has been computed but
292: before the line search routine has been applied. Allows the user to adjust the result of (usually a linear solve) that
293: determined the search direction.
295: Logically Collective
297: Input Parameters:
298: + linesearch - the `SNESLineSearch` context
299: . func - [optional] function evaluation routine, for the calling sequence see `SNESLineSearchPreCheck()`
300: - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
302: Level: intermediate
304: Note:
305: Use `SNESLineSearchSetPostCheck()` to change the step after the line search.
306: search is complete.
308: Use `SNESVISetVariableBounds()` and `SNESVISetComputeVariableBounds()` to cause `SNES` to automatically control the ranges of variables allowed.
310: .seealso: `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()`,
311: `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()
313: @*/
314: PetscErrorCode SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch, Vec, Vec, PetscBool *, void *), void *ctx)
315: {
316: PetscFunctionBegin;
318: if (func) linesearch->ops->precheck = func;
319: if (ctx) linesearch->precheckctx = ctx;
320: PetscFunctionReturn(PETSC_SUCCESS);
321: }
323: /*@C
324: SNESLineSearchGetPreCheck - Gets the pre-check function for the line search routine.
326: Input Parameter:
327: . linesearch - the `SNESLineSearch` context
329: Output Parameters:
330: + func - [optional] function evaluation routine, for calling sequence see `SNESLineSearchPreCheck`
331: - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
333: Level: intermediate
335: .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchSetPostCheck()`
336: @*/
337: PetscErrorCode SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch, Vec, Vec, PetscBool *, void *), void **ctx)
338: {
339: PetscFunctionBegin;
341: if (func) *func = linesearch->ops->precheck;
342: if (ctx) *ctx = linesearch->precheckctx;
343: PetscFunctionReturn(PETSC_SUCCESS);
344: }
346: /*@C
347: SNESLineSearchSetPostCheck - Sets a user function that is called after the line search has been applied to determine the step
348: direction and length. Allows the user a chance to change or override the decision of the line search routine
350: Logically Collective
352: Input Parameters:
353: + linesearch - the `SNESLineSearch` context
354: . func - [optional] function evaluation routine, for the calling sequence see `SNESLineSearchPostCheck`
355: - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
357: Level: intermediate
359: Notes:
360: Use `SNESLineSearchSetPreCheck()` to change the step before the line search.
361: search is complete.
363: Use `SNESVISetVariableBounds()` and `SNESVISetComputeVariableBounds()` to cause `SNES` to automatically control the ranges of variables allowed.
365: .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchGetPostCheck()`,
366: `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()
367: @*/
368: PetscErrorCode SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx)
369: {
370: PetscFunctionBegin;
372: if (func) linesearch->ops->postcheck = func;
373: if (ctx) linesearch->postcheckctx = ctx;
374: PetscFunctionReturn(PETSC_SUCCESS);
375: }
377: /*@C
378: SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine.
380: Input Parameter:
381: . linesearch - the `SNESLineSearch` context
383: Output Parameters:
384: + func - [optional] function evaluation routine, see for the calling sequence `SNESLineSearchPostCheck`
385: - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
387: Level: intermediate
389: .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`
390: @*/
391: PetscErrorCode SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx)
392: {
393: PetscFunctionBegin;
395: if (func) *func = linesearch->ops->postcheck;
396: if (ctx) *ctx = linesearch->postcheckctx;
397: PetscFunctionReturn(PETSC_SUCCESS);
398: }
400: /*@
401: SNESLineSearchPreCheck - Prepares the line search for being applied.
403: Logically Collective
405: Input Parameters:
406: + linesearch - The linesearch instance.
407: . X - The current solution
408: - Y - The step direction
410: Output Parameter:
411: . changed - Indicator that the precheck routine has changed anything
413: Level: advanced
415: Note:
416: This calls any function provided with `SNESLineSearchSetPreCheck()`
418: .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`,
419: `SNESLineSearchGetPostCheck()``
420: @*/
421: PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed)
422: {
423: PetscFunctionBegin;
424: *changed = PETSC_FALSE;
425: if (linesearch->ops->precheck) {
426: PetscUseTypeMethod(linesearch, precheck, X, Y, changed, linesearch->precheckctx);
428: }
429: PetscFunctionReturn(PETSC_SUCCESS);
430: }
432: /*@
433: SNESLineSearchPostCheck - Hook to modify step direction or updated solution after a successful linesearch
435: Logically Collective
437: Input Parameters:
438: + linesearch - The linesearch context
439: . X - The last solution
440: . Y - The step direction
441: - W - The updated solution, W = X + lambda*Y for some lambda
443: Output Parameters:
444: + changed_Y - Indicator if the direction Y has been changed.
445: - changed_W - Indicator if the new candidate solution W has been changed.
447: Level: developer
449: Note:
450: This calls any function provided with `SNESLineSearchSetPreCheck()`
452: .seealso: `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPrecheck()`, `SNESLineSearchGetPrecheck()`
453: @*/
454: PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W)
455: {
456: PetscFunctionBegin;
457: *changed_Y = PETSC_FALSE;
458: *changed_W = PETSC_FALSE;
459: if (linesearch->ops->postcheck) {
460: PetscUseTypeMethod(linesearch, postcheck, X, Y, W, changed_Y, changed_W, linesearch->postcheckctx);
463: }
464: PetscFunctionReturn(PETSC_SUCCESS);
465: }
467: /*@C
468: SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration
470: Logically Collective
472: Input Parameters:
473: + linesearch - linesearch context
474: . X - base state for this step
475: - ctx - context for this function
477: Input/Output Parameter:
478: . Y - correction, possibly modified
480: Output Parameter:
481: . changed - flag indicating that Y was modified
483: Options Database Key:
484: + -snes_linesearch_precheck_picard - activate this routine
485: - -snes_linesearch_precheck_picard_angle - angle
487: Level: advanced
489: Notes:
490: This function should be passed to `SNESLineSearchSetPreCheck()`
492: The justification for this method involves the linear convergence of a Picard iteration
493: so the Picard linearization should be provided in place of the "Jacobian". This correction
494: is generally not useful when using a Newton linearization.
496: Reference:
497: . - * - Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology.
499: .seealso: `SNESLineSearch`, `SNESSetPicard()`, `SNESGetLineSearch()`, `SNESLineSearchSetPreCheck()`
500: @*/
501: PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed, void *ctx)
502: {
503: PetscReal angle = *(PetscReal *)linesearch->precheckctx;
504: Vec Ylast;
505: PetscScalar dot;
506: PetscInt iter;
507: PetscReal ynorm, ylastnorm, theta, angle_radians;
508: SNES snes;
510: PetscFunctionBegin;
511: PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
512: PetscCall(PetscObjectQuery((PetscObject)snes, "SNESLineSearchPreCheckPicard_Ylast", (PetscObject *)&Ylast));
513: if (!Ylast) {
514: PetscCall(VecDuplicate(Y, &Ylast));
515: PetscCall(PetscObjectCompose((PetscObject)snes, "SNESLineSearchPreCheckPicard_Ylast", (PetscObject)Ylast));
516: PetscCall(PetscObjectDereference((PetscObject)Ylast));
517: }
518: PetscCall(SNESGetIterationNumber(snes, &iter));
519: if (iter < 2) {
520: PetscCall(VecCopy(Y, Ylast));
521: *changed = PETSC_FALSE;
522: PetscFunctionReturn(PETSC_SUCCESS);
523: }
525: PetscCall(VecDot(Y, Ylast, &dot));
526: PetscCall(VecNorm(Y, NORM_2, &ynorm));
527: PetscCall(VecNorm(Ylast, NORM_2, &ylastnorm));
528: /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
529: theta = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm), -1.0, 1.0));
530: angle_radians = angle * PETSC_PI / 180.;
531: if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
532: /* Modify the step Y */
533: PetscReal alpha, ydiffnorm;
534: PetscCall(VecAXPY(Ylast, -1.0, Y));
535: PetscCall(VecNorm(Ylast, NORM_2, &ydiffnorm));
536: alpha = (ydiffnorm > .001 * ylastnorm) ? ylastnorm / ydiffnorm : 1000.0;
537: PetscCall(VecCopy(Y, Ylast));
538: PetscCall(VecScale(Y, alpha));
539: 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));
540: *changed = PETSC_TRUE;
541: } else {
542: PetscCall(PetscInfo(snes, "Angle %14.12e degrees exceeds threshold %14.12e, no correction applied\n", (double)(theta * 180. / PETSC_PI), (double)angle));
543: PetscCall(VecCopy(Y, Ylast));
544: *changed = PETSC_FALSE;
545: }
546: PetscFunctionReturn(PETSC_SUCCESS);
547: }
549: /*@
550: SNESLineSearchApply - Computes the line-search update.
552: Collective
554: Input Parameters:
555: + linesearch - The linesearch context
556: - Y - The search direction
558: Input/Output Parameters:
559: + X - The current solution, on output the new solution
560: . F - The current function, on output the new function
561: - fnorm - The current norm, on output the new function norm
563: Options Database Keys:
564: + -snes_linesearch_type - basic (or equivalently none), bt, l2, cp, nleqerr, shell
565: . -snes_linesearch_monitor [:filename] - Print progress of line searches
566: . -snes_linesearch_damping - The linesearch damping parameter, default is 1.0 (no damping)
567: . -snes_linesearch_norms - Turn on/off the linesearch norms computation (SNESLineSearchSetComputeNorms())
568: . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess
569: - -snes_linesearch_max_it - The number of iterations for iterative line searches
571: Level: Intermediate
573: Notes:
574: This is typically called from within a `SNESSolve()` implementation in order to
575: help with convergence of the nonlinear method. Various `SNES` types use line searches
576: in different ways, but the overarching theme is that a line search is used to determine
577: an optimal damping parameter of a step at each iteration of the method. Each
578: application of the line search may invoke `SNESComputeFunction()` several times, and
579: therefore may be fairly expensive.
581: .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchPreCheck()`, `SNESLineSearchPostCheck()`, `SNESSolve()`, `SNESComputeFunction()`, `SNESLineSearchSetComputeNorms()`,
582: `SNESLineSearchType`, `SNESLineSearchSetType()`
583: @*/
584: PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal *fnorm, Vec Y)
585: {
586: PetscFunctionBegin;
592: linesearch->result = SNES_LINESEARCH_SUCCEEDED;
594: linesearch->vec_sol = X;
595: linesearch->vec_update = Y;
596: linesearch->vec_func = F;
598: PetscCall(SNESLineSearchSetUp(linesearch));
600: if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */
602: if (fnorm) linesearch->fnorm = *fnorm;
603: else PetscCall(VecNorm(F, NORM_2, &linesearch->fnorm));
605: PetscCall(PetscLogEventBegin(SNESLINESEARCH_Apply, linesearch, X, F, Y));
607: PetscUseTypeMethod(linesearch, apply);
609: PetscCall(PetscLogEventEnd(SNESLINESEARCH_Apply, linesearch, X, F, Y));
611: if (fnorm) *fnorm = linesearch->fnorm;
612: PetscFunctionReturn(PETSC_SUCCESS);
613: }
615: /*@
616: SNESLineSearchDestroy - Destroys the line search instance.
618: Collective
620: Input Parameter:
621: . linesearch - The linesearch context
623: Level: developer
625: .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchReset()`, `SNESDestroy()`
626: @*/
627: PetscErrorCode SNESLineSearchDestroy(SNESLineSearch *linesearch)
628: {
629: PetscFunctionBegin;
630: if (!*linesearch) PetscFunctionReturn(PETSC_SUCCESS);
632: if (--((PetscObject)(*linesearch))->refct > 0) {
633: *linesearch = NULL;
634: PetscFunctionReturn(PETSC_SUCCESS);
635: }
636: PetscCall(PetscObjectSAWsViewOff((PetscObject)*linesearch));
637: PetscCall(SNESLineSearchReset(*linesearch));
638: PetscTryTypeMethod(*linesearch, destroy);
639: PetscCall(PetscViewerDestroy(&(*linesearch)->monitor));
640: PetscCall(SNESLineSearchMonitorCancel((*linesearch)));
641: PetscCall(PetscHeaderDestroy(linesearch));
642: PetscFunctionReturn(PETSC_SUCCESS);
643: }
645: /*@
646: SNESLineSearchSetDefaultMonitor - Turns on/off printing useful information and debugging output about the line search.
648: Logically Collective
650: Input Parameters:
651: + linesearch - the linesearch object
652: - viewer - an `PETSCVIEWERASCII` `PetscViewer` or `NULL` to turn off monitor
654: Options Database Key:
655: . -snes_linesearch_monitor [:filename] - enables the monitor
657: Level: intermediate
659: Developer Note:
660: This monitor is implemented differently than the other line search monitors that are set with
661: `SNESLineSearchMonitorSet()` since it is called in many locations of the line search routines to display aspects of the
662: line search that are not visible to the other monitors.
664: .seealso: `SNESLineSearch`, `PETSCVIEWERASCII`, `SNESGetLineSearch()`, `SNESLineSearchGetDefaultMonitor()`, `PetscViewer`, `SNESLineSearchSetMonitor()`,
665: `SNESLineSearchMonitorSetFromOptions()`
666: @*/
667: PetscErrorCode SNESLineSearchSetDefaultMonitor(SNESLineSearch linesearch, PetscViewer viewer)
668: {
669: PetscFunctionBegin;
670: if (viewer) PetscCall(PetscObjectReference((PetscObject)viewer));
671: PetscCall(PetscViewerDestroy(&linesearch->monitor));
672: linesearch->monitor = viewer;
673: PetscFunctionReturn(PETSC_SUCCESS);
674: }
676: /*@
677: SNESLineSearchGetDefaultMonitor - Gets the `PetscViewer` instance for the line search monitor.
679: Logically Collective
681: Input Parameter:
682: . linesearch - linesearch context
684: Output Parameter:
685: . monitor - monitor context
687: Level: intermediate
689: .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchSetDefaultMonitor()`, `PetscViewer`
690: @*/
691: PetscErrorCode SNESLineSearchGetDefaultMonitor(SNESLineSearch linesearch, PetscViewer *monitor)
692: {
693: PetscFunctionBegin;
695: *monitor = linesearch->monitor;
696: PetscFunctionReturn(PETSC_SUCCESS);
697: }
699: /*@C
700: SNESLineSearchMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
702: Collective
704: Input Parameters:
705: + ls - `SNESLineSearch` object you wish to monitor
706: . name - the monitor type one is seeking
707: . help - message indicating what monitoring is done
708: . manual - manual page for the monitor
709: . monitor - the monitor function
710: - monitorsetup - a function that is called once ONLY if the user selected this monitor that may set additional features of the `SNESLineSearch` or `PetscViewer`
712: Level: developer
714: .seealso: `SNESLineSearch`, `SNESLineSearchSetMonitor()`, `PetscOptionsGetViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
715: `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
716: `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
717: `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
718: `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
719: `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
720: `PetscOptionsFList()`, `PetscOptionsEList()`
721: @*/
722: PetscErrorCode SNESLineSearchMonitorSetFromOptions(SNESLineSearch ls, const char name[], const char help[], const char manual[], PetscErrorCode (*monitor)(SNESLineSearch, PetscViewerAndFormat *), PetscErrorCode (*monitorsetup)(SNESLineSearch, PetscViewerAndFormat *))
723: {
724: PetscViewer viewer;
725: PetscViewerFormat format;
726: PetscBool flg;
728: PetscFunctionBegin;
729: PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)ls), ((PetscObject)ls)->options, ((PetscObject)ls)->prefix, name, &viewer, &format, &flg));
730: if (flg) {
731: PetscViewerAndFormat *vf;
732: PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf));
733: PetscCall(PetscObjectDereference((PetscObject)viewer));
734: if (monitorsetup) PetscCall((*monitorsetup)(ls, vf));
735: PetscCall(SNESLineSearchMonitorSet(ls, (PetscErrorCode(*)(SNESLineSearch, void *))monitor, vf, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy));
736: }
737: PetscFunctionReturn(PETSC_SUCCESS);
738: }
740: /*@
741: SNESLineSearchSetFromOptions - Sets options for the line search
743: Logically Collective
745: Input Parameter:
746: . linesearch - linesearch context
748: Options Database Keys:
749: + -snes_linesearch_type <type> - basic (or equivalently none), bt, l2, cp, nleqerr, shell
750: . -snes_linesearch_order <order> - 1, 2, 3. Most types only support certain orders (bt supports 2 or 3)
751: . -snes_linesearch_norms - Turn on/off the linesearch norms for the basic linesearch typem (`SNESLineSearchSetComputeNorms()`)
752: . -snes_linesearch_minlambda - The minimum step length
753: . -snes_linesearch_maxstep - The maximum step size
754: . -snes_linesearch_rtol - Relative tolerance for iterative line searches
755: . -snes_linesearch_atol - Absolute tolerance for iterative line searches
756: . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches
757: . -snes_linesearch_max_it - The number of iterations for iterative line searches
758: . -snes_linesearch_monitor [:filename] - Print progress of line searches
759: . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine
760: . -snes_linesearch_damping - The linesearch damping parameter
761: . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess.
762: . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method
763: - -snes_linesearch_precheck_picard_angle - Angle used in Picard precheck method
765: Level: intermediate
767: .seealso: `SNESLineSearch`, `SNESLineSearchCreate()`, `SNESLineSearchSetOrder()`, `SNESLineSearchSetType()`, `SNESLineSearchSetTolerances()`, `SNESLineSearchSetDamping()`, `SNESLineSearchPreCheckPicard()`,
768: `SNESLineSearchType`, `SNESLineSearchSetComputeNorms()`
769: @*/
770: PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch)
771: {
772: const char *deft = SNESLINESEARCHBASIC;
773: char type[256];
774: PetscBool flg, set;
775: PetscViewer viewer;
777: PetscFunctionBegin;
778: PetscCall(SNESLineSearchRegisterAll());
780: PetscObjectOptionsBegin((PetscObject)linesearch);
781: if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name;
782: PetscCall(PetscOptionsFList("-snes_linesearch_type", "Linesearch type", "SNESLineSearchSetType", SNESLineSearchList, deft, type, 256, &flg));
783: if (flg) {
784: PetscCall(SNESLineSearchSetType(linesearch, type));
785: } else if (!((PetscObject)linesearch)->type_name) {
786: PetscCall(SNESLineSearchSetType(linesearch, deft));
787: }
789: PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)linesearch), ((PetscObject)linesearch)->options, ((PetscObject)linesearch)->prefix, "-snes_linesearch_monitor", &viewer, NULL, &set));
790: if (set) {
791: PetscCall(SNESLineSearchSetDefaultMonitor(linesearch, viewer));
792: PetscCall(PetscViewerDestroy(&viewer));
793: }
794: PetscCall(SNESLineSearchMonitorSetFromOptions(linesearch, "-snes_linesearch_monitor_solution_update", "View correction at each iteration", "SNESLineSearchMonitorSolutionUpdate", SNESLineSearchMonitorSolutionUpdate, NULL));
796: /* tolerances */
797: PetscCall(PetscOptionsReal("-snes_linesearch_minlambda", "Minimum step length", "SNESLineSearchSetTolerances", linesearch->steptol, &linesearch->steptol, NULL));
798: PetscCall(PetscOptionsReal("-snes_linesearch_maxstep", "Maximum step size", "SNESLineSearchSetTolerances", linesearch->maxstep, &linesearch->maxstep, NULL));
799: PetscCall(PetscOptionsReal("-snes_linesearch_rtol", "Relative tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->rtol, &linesearch->rtol, NULL));
800: PetscCall(PetscOptionsReal("-snes_linesearch_atol", "Absolute tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->atol, &linesearch->atol, NULL));
801: PetscCall(PetscOptionsReal("-snes_linesearch_ltol", "Change in lambda tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->ltol, &linesearch->ltol, NULL));
802: PetscCall(PetscOptionsInt("-snes_linesearch_max_it", "Maximum iterations for iterative line searches", "SNESLineSearchSetTolerances", linesearch->max_its, &linesearch->max_its, NULL));
804: /* damping parameters */
805: PetscCall(PetscOptionsReal("-snes_linesearch_damping", "Line search damping and initial step guess", "SNESLineSearchSetDamping", linesearch->damping, &linesearch->damping, NULL));
807: PetscCall(PetscOptionsBool("-snes_linesearch_keeplambda", "Use previous lambda as damping", "SNESLineSearchSetKeepLambda", linesearch->keeplambda, &linesearch->keeplambda, NULL));
809: /* precheck */
810: PetscCall(PetscOptionsBool("-snes_linesearch_precheck_picard", "Use a correction that sometimes improves convergence of Picard iteration", "SNESLineSearchPreCheckPicard", flg, &flg, &set));
811: if (set) {
812: if (flg) {
813: linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */
815: 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));
816: PetscCall(SNESLineSearchSetPreCheck(linesearch, SNESLineSearchPreCheckPicard, &linesearch->precheck_picard_angle));
817: } else {
818: PetscCall(SNESLineSearchSetPreCheck(linesearch, NULL, NULL));
819: }
820: }
821: PetscCall(PetscOptionsInt("-snes_linesearch_order", "Order of approximation used in the line search", "SNESLineSearchSetOrder", linesearch->order, &linesearch->order, NULL));
822: PetscCall(PetscOptionsBool("-snes_linesearch_norms", "Compute final norms in line search", "SNESLineSearchSetComputeNorms", linesearch->norms, &linesearch->norms, NULL));
824: PetscTryTypeMethod(linesearch, setfromoptions, PetscOptionsObject);
826: PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)linesearch, PetscOptionsObject));
827: PetscOptionsEnd();
828: PetscFunctionReturn(PETSC_SUCCESS);
829: }
831: /*@
832: SNESLineSearchView - Prints useful information about the line search
834: Logically Collective
836: Input Parameters:
837: + linesearch - linesearch context
838: - viewer - the viewer to display the line search information
840: Level: intermediate
842: .seealso: `SNESLineSearch`, `PetscViewer`, `SNESLineSearchCreate()`
843: @*/
844: PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer)
845: {
846: PetscBool iascii;
848: PetscFunctionBegin;
850: if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch), &viewer));
852: PetscCheckSameComm(linesearch, 1, viewer, 2);
854: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
855: if (iascii) {
856: PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)linesearch, viewer));
857: PetscCall(PetscViewerASCIIPushTab(viewer));
858: PetscTryTypeMethod(linesearch, view, viewer);
859: PetscCall(PetscViewerASCIIPopTab(viewer));
860: PetscCall(PetscViewerASCIIPrintf(viewer, " maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep, (double)linesearch->steptol));
861: PetscCall(PetscViewerASCIIPrintf(viewer, " tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol, (double)linesearch->atol, (double)linesearch->ltol));
862: PetscCall(PetscViewerASCIIPrintf(viewer, " maximum iterations=%" PetscInt_FMT "\n", linesearch->max_its));
863: if (linesearch->ops->precheck) {
864: if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) {
865: PetscCall(PetscViewerASCIIPrintf(viewer, " using precheck step to speed up Picard convergence\n"));
866: } else {
867: PetscCall(PetscViewerASCIIPrintf(viewer, " using user-defined precheck step\n"));
868: }
869: }
870: if (linesearch->ops->postcheck) PetscCall(PetscViewerASCIIPrintf(viewer, " using user-defined postcheck step\n"));
871: }
872: PetscFunctionReturn(PETSC_SUCCESS);
873: }
875: /*@C
876: SNESLineSearchGetType - Gets the linesearch type
878: Logically Collective
880: Input Parameter:
881: . linesearch - linesearch context
883: Output Parameter:
884: . type - The type of line search, or `NULL` if not set
886: Level: intermediate
888: .seealso: `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchType`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchSetType()`
889: @*/
890: PetscErrorCode SNESLineSearchGetType(SNESLineSearch linesearch, SNESLineSearchType *type)
891: {
892: PetscFunctionBegin;
895: *type = ((PetscObject)linesearch)->type_name;
896: PetscFunctionReturn(PETSC_SUCCESS);
897: }
899: /*@C
900: SNESLineSearchSetType - Sets the linesearch type
902: Logically Collective
904: Input Parameters:
905: + linesearch - linesearch context
906: - type - The type of line search to be used
908: Available Types:
909: + `SNESLINESEARCHBASIC` - (or equivalently `SNESLINESEARCHNONE`) Simple damping line search, defaults to using the full Newton step
910: . `SNESLINESEARCHBT` - Backtracking line search over the L2 norm of the function
911: . `SNESLINESEARCHL2` - Secant line search over the L2 norm of the function
912: . `SNESLINESEARCHCP` - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x)
913: . `SNESLINESEARCHNLEQERR` - Affine-covariant error-oriented linesearch
914: - `SNESLINESEARCHSHELL` - User provided `SNESLineSearch` implementation
916: Options Database Key:
917: . -snes_linesearch_type <type> - basic (or equivalently none), bt, l2, cp, nleqerr, shell
919: Level: intermediate
921: .seealso: `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchType`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchGetType()`
922: @*/
923: PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type)
924: {
925: PetscBool match;
926: PetscErrorCode (*r)(SNESLineSearch);
928: PetscFunctionBegin;
932: PetscCall(PetscObjectTypeCompare((PetscObject)linesearch, type, &match));
933: if (match) PetscFunctionReturn(PETSC_SUCCESS);
935: PetscCall(PetscFunctionListFind(SNESLineSearchList, type, &r));
936: PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested Line Search type %s", type);
937: /* Destroy the previous private linesearch context */
938: if (linesearch->ops->destroy) {
939: PetscCall((*(linesearch)->ops->destroy)(linesearch));
940: linesearch->ops->destroy = NULL;
941: }
942: /* Reinitialize function pointers in SNESLineSearchOps structure */
943: linesearch->ops->apply = NULL;
944: linesearch->ops->view = NULL;
945: linesearch->ops->setfromoptions = NULL;
946: linesearch->ops->destroy = NULL;
948: PetscCall(PetscObjectChangeTypeName((PetscObject)linesearch, type));
949: PetscCall((*r)(linesearch));
950: PetscFunctionReturn(PETSC_SUCCESS);
951: }
953: /*@
954: SNESLineSearchSetSNES - Sets the `SNES` for the linesearch for function evaluation.
956: Input Parameters:
957: + linesearch - linesearch context
958: - snes - The snes instance
960: Level: developer
962: Note:
963: This happens automatically when the line search is obtained/created with
964: `SNESGetLineSearch()`. This routine is therefore mainly called within `SNES`
965: implementations.
967: .seealso: `SNESLineSearch`, `SNESLineSearchGetSNES()`, `SNESLineSearchSetVecs()`, `SNES`
968: @*/
969: PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes)
970: {
971: PetscFunctionBegin;
974: linesearch->snes = snes;
975: PetscFunctionReturn(PETSC_SUCCESS);
976: }
978: /*@
979: SNESLineSearchGetSNES - Gets the `SNES` instance associated with the line search.
980: Having an associated `SNES` is necessary because most line search implementations must be able to
981: evaluate the function using `SNESComputeFunction()` for the associated `SNES`. This routine
982: is used in the line search implementations when one must get this associated `SNES` instance.
984: Not Collective
986: Input Parameter:
987: . linesearch - linesearch context
989: Output Parameter:
990: . snes - The `SNES` instance
992: Level: developer
994: .seealso: `SNESLineSearch`, `SNESType`, `SNESLineSearchGetSNES()`, `SNESLineSearchSetVecs()`, `SNES`
995: @*/
996: PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes)
997: {
998: PetscFunctionBegin;
1001: *snes = linesearch->snes;
1002: PetscFunctionReturn(PETSC_SUCCESS);
1003: }
1005: /*@
1006: SNESLineSearchGetLambda - Gets the last linesearch steplength discovered.
1008: Not Collective
1010: Input Parameter:
1011: . linesearch - linesearch context
1013: Output Parameter:
1014: . lambda - The last steplength computed during `SNESLineSearchApply()`
1016: Level: advanced
1018: Note:
1019: This is useful in methods where the solver is ill-scaled and
1020: requires some adaptive notion of the difference in scale between the
1021: solution and the function. For instance, `SNESQN` may be scaled by the
1022: line search lambda using the argument -snes_qn_scaling ls.
1024: .seealso: `SNESLineSearch`, `SNESLineSearchSetLambda()`, `SNESLineSearchGetDamping()`, `SNESLineSearchApply()`
1025: @*/
1026: PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch linesearch, PetscReal *lambda)
1027: {
1028: PetscFunctionBegin;
1031: *lambda = linesearch->lambda;
1032: PetscFunctionReturn(PETSC_SUCCESS);
1033: }
1035: /*@
1036: SNESLineSearchSetLambda - Sets the linesearch steplength
1038: Input Parameters:
1039: + linesearch - linesearch context
1040: - lambda - The last steplength.
1042: Level: advanced
1044: Note:
1045: This routine is typically used within implementations of `SNESLineSearchApply()`
1046: to set the final steplength. This routine (and `SNESLineSearchGetLambda()`) were
1047: added in order to facilitate Quasi-Newton methods that use the previous steplength
1048: as an inner scaling parameter.
1050: .seealso: `SNESLineSearch`, `SNESLineSearchGetLambda()`
1051: @*/
1052: PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
1053: {
1054: PetscFunctionBegin;
1056: linesearch->lambda = lambda;
1057: PetscFunctionReturn(PETSC_SUCCESS);
1058: }
1060: /*@
1061: SNESLineSearchGetTolerances - Gets the tolerances for the linesearch. These include
1062: tolerances for the relative and absolute change in the function norm, the change
1063: in lambda for iterative line searches, the minimum steplength, the maximum steplength,
1064: and the maximum number of iterations the line search procedure may take.
1066: Not Collective
1068: Input Parameter:
1069: . linesearch - linesearch context
1071: Output Parameters:
1072: + steptol - The minimum steplength
1073: . maxstep - The maximum steplength
1074: . rtol - The relative tolerance for iterative line searches
1075: . atol - The absolute tolerance for iterative line searches
1076: . ltol - The change in lambda tolerance for iterative line searches
1077: - max_it - The maximum number of iterations of the line search
1079: Level: intermediate
1081: Note:
1082: Different line searches may implement these parameters slightly differently as
1083: the type requires.
1085: .seealso: `SNESLineSearch`, `SNESLineSearchSetTolerances()`
1086: @*/
1087: PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch linesearch, PetscReal *steptol, PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
1088: {
1089: PetscFunctionBegin;
1091: if (steptol) {
1093: *steptol = linesearch->steptol;
1094: }
1095: if (maxstep) {
1097: *maxstep = linesearch->maxstep;
1098: }
1099: if (rtol) {
1101: *rtol = linesearch->rtol;
1102: }
1103: if (atol) {
1105: *atol = linesearch->atol;
1106: }
1107: if (ltol) {
1109: *ltol = linesearch->ltol;
1110: }
1111: if (max_its) {
1113: *max_its = linesearch->max_its;
1114: }
1115: PetscFunctionReturn(PETSC_SUCCESS);
1116: }
1118: /*@
1119: SNESLineSearchSetTolerances - Gets the tolerances for the linesearch. These include
1120: tolerances for the relative and absolute change in the function norm, the change
1121: in lambda for iterative line searches, the minimum steplength, the maximum steplength,
1122: and the maximum number of iterations the line search procedure may take.
1124: Collective
1126: Input Parameters:
1127: + linesearch - linesearch context
1128: . steptol - The minimum steplength
1129: . maxstep - The maximum steplength
1130: . rtol - The relative tolerance for iterative line searches
1131: . atol - The absolute tolerance for iterative line searches
1132: . ltol - The change in lambda tolerance for iterative line searches
1133: - max_it - The maximum number of iterations of the line search
1135: Level: intermediate
1137: Note:
1138: The user may choose to not set any of the tolerances using `PETSC_DEFAULT` in
1139: place of an argument.
1141: .seealso: `SNESLineSearch`, `SNESLineSearchGetTolerances()`
1142: @*/
1143: PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch linesearch, PetscReal steptol, PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its)
1144: {
1145: PetscFunctionBegin;
1154: if (steptol != (PetscReal)PETSC_DEFAULT) {
1155: PetscCheck(steptol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Minimum step length %14.12e must be non-negative", (double)steptol);
1156: linesearch->steptol = steptol;
1157: }
1159: if (maxstep != (PetscReal)PETSC_DEFAULT) {
1160: PetscCheck(maxstep >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum step length %14.12e must be non-negative", (double)maxstep);
1161: linesearch->maxstep = maxstep;
1162: }
1164: if (rtol != (PetscReal)PETSC_DEFAULT) {
1165: 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);
1166: linesearch->rtol = rtol;
1167: }
1169: if (atol != (PetscReal)PETSC_DEFAULT) {
1170: PetscCheck(atol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Absolute tolerance %14.12e must be non-negative", (double)atol);
1171: linesearch->atol = atol;
1172: }
1174: if (ltol != (PetscReal)PETSC_DEFAULT) {
1175: PetscCheck(ltol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Lambda tolerance %14.12e must be non-negative", (double)ltol);
1176: linesearch->ltol = ltol;
1177: }
1179: if (max_its != PETSC_DEFAULT) {
1180: PetscCheck(max_its >= 0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of iterations %" PetscInt_FMT " must be non-negative", max_its);
1181: linesearch->max_its = max_its;
1182: }
1183: PetscFunctionReturn(PETSC_SUCCESS);
1184: }
1186: /*@
1187: SNESLineSearchGetDamping - Gets the line search damping parameter.
1189: Input Parameter:
1190: . linesearch - linesearch context
1192: Output Parameter:
1193: . damping - The damping parameter
1195: Level: advanced
1197: .seealso: `SNESLineSearchGetStepTolerance()`, `SNESQN`
1198: @*/
1200: PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch linesearch, PetscReal *damping)
1201: {
1202: PetscFunctionBegin;
1205: *damping = linesearch->damping;
1206: PetscFunctionReturn(PETSC_SUCCESS);
1207: }
1209: /*@
1210: SNESLineSearchSetDamping - Sets the line search damping parameter.
1212: Input Parameters:
1213: + linesearch - linesearch context
1214: - damping - The damping parameter
1216: Options Database Key:
1217: . -snes_linesearch_damping <damping> - the damping value
1219: Level: intermediate
1221: Note:
1222: The `SNESLINESEARCHNONE` line search merely takes the update step scaled by the damping parameter.
1223: The use of the damping parameter in the l2 and cp line searches is much more subtle;
1224: it is used as a starting point in calculating the secant step. However, the eventual
1225: step may be of greater length than the damping parameter. In the bt line search it is
1226: used as the maximum possible step length, as the bt line search only backtracks.
1228: .seealso: `SNESLineSearch`, `SNESLineSearchGetDamping()`
1229: @*/
1230: PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch linesearch, PetscReal damping)
1231: {
1232: PetscFunctionBegin;
1234: linesearch->damping = damping;
1235: PetscFunctionReturn(PETSC_SUCCESS);
1236: }
1238: /*@
1239: SNESLineSearchGetOrder - Gets the line search approximation order.
1241: Input Parameter:
1242: . linesearch - linesearch context
1244: Output Parameter:
1245: . order - The order
1247: Possible Values for order:
1248: + 1 or `SNES_LINESEARCH_ORDER_LINEAR` - linear order
1249: . 2 or `SNES_LINESEARCH_ORDER_QUADRATIC` - quadratic order
1250: - 3 or `SNES_LINESEARCH_ORDER_CUBIC` - cubic order
1252: Level: intermediate
1254: .seealso: `SNESLineSearch`, `SNESLineSearchSetOrder()`
1255: @*/
1257: PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch linesearch, PetscInt *order)
1258: {
1259: PetscFunctionBegin;
1262: *order = linesearch->order;
1263: PetscFunctionReturn(PETSC_SUCCESS);
1264: }
1266: /*@
1267: SNESLineSearchSetOrder - Sets the maximum order of the polynomial fit used in the line search
1269: Input Parameters:
1270: + linesearch - linesearch context
1271: - order - The damping parameter
1273: Level: intermediate
1275: Possible Values for order:
1276: + 1 or `SNES_LINESEARCH_ORDER_LINEAR` - linear order
1277: . 2 or `SNES_LINESEARCH_ORDER_QUADRATIC` - quadratic order
1278: - 3 or `SNES_LINESEARCH_ORDER_CUBIC` - cubic order
1280: Notes:
1281: Variable orders are supported by the following line searches:
1282: + bt - cubic and quadratic
1283: - cp - linear and quadratic
1285: .seealso: `SNESLineSearch`, `SNESLineSearchGetOrder()`, `SNESLineSearchSetDamping()`
1286: @*/
1287: PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch linesearch, PetscInt order)
1288: {
1289: PetscFunctionBegin;
1291: linesearch->order = order;
1292: PetscFunctionReturn(PETSC_SUCCESS);
1293: }
1295: /*@
1296: SNESLineSearchGetNorms - Gets the norms for for X, Y, and F.
1298: Not Collective
1300: Input Parameter:
1301: . linesearch - linesearch context
1303: Output Parameters:
1304: + xnorm - The norm of the current solution
1305: . fnorm - The norm of the current function
1306: - ynorm - The norm of the current update
1308: Level: developer
1310: .seealso: `SNESLineSearch`, `SNESLineSearchSetNorms()` `SNESLineSearchGetVecs()`
1311: @*/
1312: PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal *xnorm, PetscReal *fnorm, PetscReal *ynorm)
1313: {
1314: PetscFunctionBegin;
1316: if (xnorm) *xnorm = linesearch->xnorm;
1317: if (fnorm) *fnorm = linesearch->fnorm;
1318: if (ynorm) *ynorm = linesearch->ynorm;
1319: PetscFunctionReturn(PETSC_SUCCESS);
1320: }
1322: /*@
1323: SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F.
1325: Collective
1327: Input Parameters:
1328: + linesearch - linesearch context
1329: . xnorm - The norm of the current solution
1330: . fnorm - The norm of the current function
1331: - ynorm - The norm of the current update
1333: Level: developer
1335: .seealso: `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()`
1336: @*/
1337: PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
1338: {
1339: PetscFunctionBegin;
1341: linesearch->xnorm = xnorm;
1342: linesearch->fnorm = fnorm;
1343: linesearch->ynorm = ynorm;
1344: PetscFunctionReturn(PETSC_SUCCESS);
1345: }
1347: /*@
1348: SNESLineSearchComputeNorms - Computes the norms of X, F, and Y.
1350: Input Parameter:
1351: . linesearch - linesearch context
1353: Options Database Key:
1354: . -snes_linesearch_norms - turn norm computation on or off
1356: Level: intermediate
1358: .seealso: `SNESLineSearch`, `SNESLineSearchGetNorms`, `SNESLineSearchSetNorms()`, `SNESLineSearchSetComputeNorms()`
1359: @*/
1360: PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
1361: {
1362: SNES snes;
1364: PetscFunctionBegin;
1365: if (linesearch->norms) {
1366: if (linesearch->ops->vinorm) {
1367: PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
1368: PetscCall(VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
1369: PetscCall(VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm));
1370: PetscCall((*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm));
1371: } else {
1372: PetscCall(VecNormBegin(linesearch->vec_func, NORM_2, &linesearch->fnorm));
1373: PetscCall(VecNormBegin(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
1374: PetscCall(VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm));
1375: PetscCall(VecNormEnd(linesearch->vec_func, NORM_2, &linesearch->fnorm));
1376: PetscCall(VecNormEnd(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
1377: PetscCall(VecNormEnd(linesearch->vec_update, NORM_2, &linesearch->ynorm));
1378: }
1379: }
1380: PetscFunctionReturn(PETSC_SUCCESS);
1381: }
1383: /*@
1384: SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search.
1386: Input Parameters:
1387: + linesearch - linesearch context
1388: - flg - indicates whether or not to compute norms
1390: Options Database Key:
1391: . -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic (none) linesearch
1393: Level: intermediate
1395: Note:
1396: This is most relevant to the `SNESLINESEARCHBASIC` (or equivalently `SNESLINESEARCHNONE`) line search type since most line searches have a stopping criteria involving the norm.
1398: .seealso: `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetNorms()`, `SNESLineSearchComputeNorms()`, `SNESLINESEARCHBASIC`
1399: @*/
1400: PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg)
1401: {
1402: PetscFunctionBegin;
1403: linesearch->norms = flg;
1404: PetscFunctionReturn(PETSC_SUCCESS);
1405: }
1407: /*@
1408: SNESLineSearchGetVecs - Gets the vectors from the `SNESLineSearch` context
1410: Not Collective on but the vectors are parallel
1412: Input Parameter:
1413: . linesearch - linesearch context
1415: Output Parameters:
1416: + X - Solution vector
1417: . F - Function vector
1418: . Y - Search direction vector
1419: . W - Solution work vector
1420: - G - Function work vector
1422: Level: advanced
1424: Notes:
1425: At the beginning of a line search application, `X` should contain a
1426: solution and the vector `F` the function computed at `X`. At the end of the
1427: line search application, `X` should contain the new solution, and `F` the
1428: function evaluated at the new solution.
1430: These vectors are owned by the `SNESLineSearch` and should not be destroyed by the caller
1432: .seealso: `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()`
1433: @*/
1434: PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch, Vec *X, Vec *F, Vec *Y, Vec *W, Vec *G)
1435: {
1436: PetscFunctionBegin;
1438: if (X) {
1440: *X = linesearch->vec_sol;
1441: }
1442: if (F) {
1444: *F = linesearch->vec_func;
1445: }
1446: if (Y) {
1448: *Y = linesearch->vec_update;
1449: }
1450: if (W) {
1452: *W = linesearch->vec_sol_new;
1453: }
1454: if (G) {
1456: *G = linesearch->vec_func_new;
1457: }
1458: PetscFunctionReturn(PETSC_SUCCESS);
1459: }
1461: /*@
1462: SNESLineSearchSetVecs - Sets the vectors on the `SNESLineSearch` context
1464: Logically Collective
1466: Input Parameters:
1467: + linesearch - linesearch context
1468: . X - Solution vector
1469: . F - Function vector
1470: . Y - Search direction vector
1471: . W - Solution work vector
1472: - G - Function work vector
1474: Level: advanced
1476: .seealso: `SNESLineSearch`, `SNESLineSearchSetNorms()`, `SNESLineSearchGetVecs()`
1477: @*/
1478: PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch, Vec X, Vec F, Vec Y, Vec W, Vec G)
1479: {
1480: PetscFunctionBegin;
1482: if (X) {
1484: linesearch->vec_sol = X;
1485: }
1486: if (F) {
1488: linesearch->vec_func = F;
1489: }
1490: if (Y) {
1492: linesearch->vec_update = Y;
1493: }
1494: if (W) {
1496: linesearch->vec_sol_new = W;
1497: }
1498: if (G) {
1500: linesearch->vec_func_new = G;
1501: }
1502: PetscFunctionReturn(PETSC_SUCCESS);
1503: }
1505: /*@C
1506: SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1507: `SNESLineSearch` options in the database.
1509: Logically Collective
1511: Input Parameters:
1512: + linesearch - the `SNESLineSearch` context
1513: - prefix - the prefix to prepend to all option names
1515: Level: advanced
1517: Note:
1518: A hyphen (-) must NOT be given at the beginning of the prefix name.
1519: The first character of all runtime options is AUTOMATICALLY the hyphen.
1521: .seealso: `SNESLineSearch()`, `SNESLineSearchSetFromOptions()`, `SNESGetOptionsPrefix()`
1522: @*/
1523: PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch, const char prefix[])
1524: {
1525: PetscFunctionBegin;
1527: PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)linesearch, prefix));
1528: PetscFunctionReturn(PETSC_SUCCESS);
1529: }
1531: /*@C
1532: SNESLineSearchGetOptionsPrefix - Gets the prefix used for searching for all
1533: SNESLineSearch options in the database.
1535: Not Collective
1537: Input Parameter:
1538: . linesearch - the `SNESLineSearch` context
1540: Output Parameter:
1541: . prefix - pointer to the prefix string used
1543: Level: advanced
1545: Fortran Note:
1546: The user should pass in a string 'prefix' of
1547: sufficient length to hold the prefix.
1549: .seealso: `SNESLineSearch`, `SNESAppendOptionsPrefix()`
1550: @*/
1551: PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch, const char *prefix[])
1552: {
1553: PetscFunctionBegin;
1555: PetscCall(PetscObjectGetOptionsPrefix((PetscObject)linesearch, prefix));
1556: PetscFunctionReturn(PETSC_SUCCESS);
1557: }
1559: /*@C
1560: SNESLineSearchSetWorkVecs - Sets work vectors for the line search.
1562: Input Parameters:
1563: + linesearch - the `SNESLineSearch` context
1564: - nwork - the number of work vectors
1566: Level: developer
1568: .seealso: `SNESLineSearch`, `SNESSetWorkVecs()`
1569: @*/
1570: PetscErrorCode SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork)
1571: {
1572: PetscFunctionBegin;
1573: PetscCheck(linesearch->vec_sol, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
1574: PetscCall(VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work));
1575: PetscFunctionReturn(PETSC_SUCCESS);
1576: }
1578: /*@
1579: SNESLineSearchGetReason - Gets the success/failure status of the last line search application
1581: Input Parameter:
1582: . linesearch - linesearch context
1584: Output Parameter:
1585: . result - The success or failure status
1587: Level: developer
1589: Note:
1590: This is typically called after `SNESLineSearchApply()` in order to determine if the line-search failed
1591: (and set the SNES convergence accordingly).
1593: .seealso: `SNESLineSearch`, `SNESLineSearchSetReason()`, `SNESLineSearchReason`
1594: @*/
1595: PetscErrorCode SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result)
1596: {
1597: PetscFunctionBegin;
1600: *result = linesearch->result;
1601: PetscFunctionReturn(PETSC_SUCCESS);
1602: }
1604: /*@
1605: SNESLineSearchSetReason - Sets the success/failure status of the last line search application
1607: Input Parameters:
1608: + linesearch - linesearch context
1609: - result - The success or failure status
1611: Level: developer
1613: Note:
1614: This is typically called in a `SNESLineSearchApply()` or a `SNESLINESEARCHSHELL` implementation to set
1615: the success or failure of the line search method.
1617: .seealso: `SNESLineSearch`, `SNESLineSearchReason`, `SNESLineSearchGetSResult()`
1618: @*/
1619: PetscErrorCode SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result)
1620: {
1621: PetscFunctionBegin;
1623: linesearch->result = result;
1624: PetscFunctionReturn(PETSC_SUCCESS);
1625: }
1627: /*@C
1628: SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.
1630: Logically Collective
1632: Input Parameters:
1633: + snes - nonlinear context obtained from `SNESCreate()`
1634: . projectfunc - function for projecting the function to the bounds
1635: - normfunc - function for computing the norm of an active set
1637: Calling sequence of `projectfunc`:
1638: .vb
1639: PetscErrorCode projectfunc(SNES snes, Vec X)
1640: .ve
1641: + snes - nonlinear context
1642: - X - current solution, store the projected solution here
1644: Calling sequence of `normfunc`:
1645: .vb
1646: PetscErrorCode normfunc(SNES snes, Vec X, Vec F, PetscScalar *fnorm)
1647: .ve
1648: + snes - nonlinear context
1649: . X - current solution
1650: . F - current residual
1651: - fnorm - VI-specific norm of the function
1653: Level: advanced
1655: Note:
1656: The VI solvers require projection of the solution to the feasible set. `projectfunc` should implement this.
1658: The VI solvers require special evaluation of the function norm such that the norm is only calculated
1659: on the inactive set. This should be implemented by `normfunc`.
1661: .seealso: `SNESLineSearch`, `SNESLineSearchGetVIFunctions()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()`
1662: @*/
1663: PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc)
1664: {
1665: PetscFunctionBegin;
1667: if (projectfunc) linesearch->ops->viproject = projectfunc;
1668: if (normfunc) linesearch->ops->vinorm = normfunc;
1669: PetscFunctionReturn(PETSC_SUCCESS);
1670: }
1672: /*@C
1673: SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.
1675: Not Collective
1677: Input Parameter:
1678: . linesearch - the line search context, obtain with `SNESGetLineSearch()`
1680: Output Parameters:
1681: + projectfunc - function for projecting the function to the bounds
1682: - normfunc - function for computing the norm of an active set
1684: Level: advanced
1686: .seealso: `SNESLineSearch`, `SNESLineSearchSetVIFunctions()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()`
1687: @*/
1688: PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc)
1689: {
1690: PetscFunctionBegin;
1691: if (projectfunc) *projectfunc = linesearch->ops->viproject;
1692: if (normfunc) *normfunc = linesearch->ops->vinorm;
1693: PetscFunctionReturn(PETSC_SUCCESS);
1694: }
1696: /*@C
1697: SNESLineSearchRegister - register a line search method
1699: Level: advanced
1701: .seealso: `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchSetType()`
1702: @*/
1703: PetscErrorCode SNESLineSearchRegister(const char sname[], PetscErrorCode (*function)(SNESLineSearch))
1704: {
1705: PetscFunctionBegin;
1706: PetscCall(SNESInitializePackage());
1707: PetscCall(PetscFunctionListAdd(&SNESLineSearchList, sname, function));
1708: PetscFunctionReturn(PETSC_SUCCESS);
1709: }