Actual source code: linesearch.c
petsc-3.11.4 2019-09-28
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 on SNESLineSearch
14: Input Parameters:
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: Notes:
23: There is no way to clear one specific monitor from a SNESLineSearch object.
25: This does not clear the monitor set with SNESLineSearchSetDefaultMonitor() use SNESLineSearchSetDefaultMonitor(ls,NULL) to cancel
26: that one.
28: Level: intermediate
30: .keywords: SNESLineSearch, nonlinear, set, monitor
32: .seealso: SNESGetLineSearch(), SNESLineSearchMonitorDefault(), SNESLineSearchMonitorSet()
33: @*/
34: PetscErrorCode SNESLineSearchMonitorCancel(SNESLineSearch ls)
35: {
37: PetscInt i;
41: for (i=0; i<ls->numbermonitors; i++) {
42: if (ls->monitordestroy[i]) {
43: (*ls->monitordestroy[i])(&ls->monitorcontext[i]);
44: }
45: }
46: ls->numbermonitors = 0;
47: return(0);
48: }
50: /*@
51: SNESLineSearchMonitor - runs the user provided monitor routines, if they exist
53: Collective on SNES
55: Input Parameters:
56: . ls - the linesearch object
58: Notes:
59: This routine is called by the SNES implementations.
60: It does not typically need to be called by the user.
62: Level: developer
64: .seealso: SNESGetLineSearch(), SNESLineSearchMonitorSet()
65: @*/
66: PetscErrorCode SNESLineSearchMonitor(SNESLineSearch ls)
67: {
69: PetscInt i,n = ls->numbermonitors;
72: for (i=0; i<n; i++) {
73: (*ls->monitorftns[i])(ls,ls->monitorcontext[i]);
74: }
75: return(0);
76: }
78: /*@C
79: SNESLineSearchMonitorSet - Sets an ADDITIONAL function that is to be used at every
80: iteration of the nonlinear solver to display the iteration's
81: progress.
83: Logically Collective on SNESLineSearch
85: Input Parameters:
86: + ls - the SNESLineSearch context
87: . f - the monitor function
88: . mctx - [optional] user-defined context for private data for the
89: monitor routine (use NULL if no context is desired)
90: - monitordestroy - [optional] routine that frees monitor context
91: (may be NULL)
93: Notes:
94: Several different monitoring routines may be set by calling
95: SNESLineSearchMonitorSet() multiple times; all will be called in the
96: order in which they were set.
98: Fortran Notes:
99: Only a single monitor function can be set for each SNESLineSearch object
101: Level: intermediate
103: .keywords: SNESLineSearch, nonlinear, set, monitor
105: .seealso: SNESGetLineSearch(), SNESLineSearchMonitorDefault(), SNESLineSearchMonitorCancel()
106: @*/
107: PetscErrorCode SNESLineSearchMonitorSet(SNESLineSearch ls,PetscErrorCode (*f)(SNESLineSearch,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
108: {
110: PetscInt i;
111: PetscBool identical;
115: for (i=0; i<ls->numbermonitors;i++) {
116: PetscMonitorCompare((PetscErrorCode (*)(void))f,mctx,monitordestroy,(PetscErrorCode (*)(void))ls->monitorftns[i],ls->monitorcontext[i],ls->monitordestroy[i],&identical);
117: if (identical) return(0);
118: }
119: if (ls->numbermonitors >= MAXSNESLSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
120: ls->monitorftns[ls->numbermonitors] = f;
121: ls->monitordestroy[ls->numbermonitors] = monitordestroy;
122: ls->monitorcontext[ls->numbermonitors++] = (void*)mctx;
123: return(0);
124: }
126: /*@C
127: SNESLineSearchMonitorSolutionUpdate - Monitors each update a new function value the linesearch tries
129: Collective on SNESLineSearch
131: Input Parameters:
132: + ls - the SNES linesearch object
133: - vf - the context for the monitor, in this case it is an ASCII PetscViewer and format
135: Level: intermediate
137: .keywords: SNES, nonlinear, default, monitor, norm
139: .seealso: SNESGetLineSearch(), SNESMonitorSet(), SNESMonitorSolution()
140: @*/
141: PetscErrorCode SNESLineSearchMonitorSolutionUpdate(SNESLineSearch ls,PetscViewerAndFormat *vf)
142: {
144: PetscViewer viewer = vf->viewer;
145: Vec Y,W,G;
148: SNESLineSearchGetVecs(ls,NULL,NULL,&Y,&W,&G);
149: PetscViewerPushFormat(viewer,vf->format);
150: PetscViewerASCIIPrintf(viewer,"LineSearch attempted update to solution \n");
151: VecView(Y,viewer);
152: PetscViewerASCIIPrintf(viewer,"LineSearch attempted new solution \n");
153: VecView(W,viewer);
154: PetscViewerASCIIPrintf(viewer,"LineSearch attempted updated function value\n");
155: VecView(G,viewer);
156: PetscViewerPopFormat(viewer);
157: return(0);
158: }
160: /*@
161: SNESLineSearchCreate - Creates the line search context.
163: Logically Collective on Comm
165: Input Parameters:
166: . comm - MPI communicator for the line search (typically from the associated SNES context).
168: Output Parameters:
169: . outlinesearch - the new linesearch context
171: Level: developer
173: Notes:
174: The preferred calling sequence for users is to use SNESGetLineSearch() to acquire the SNESLineSearch instance
175: already associated with the SNES. This function is for developer use.
177: .keywords: LineSearch, create, context
179: .seealso: LineSearchDestroy(), SNESGetLineSearch()
180: @*/
182: PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch)
183: {
185: SNESLineSearch linesearch;
189: SNESInitializePackage();
190: *outlinesearch = NULL;
192: PetscHeaderCreate(linesearch,SNESLINESEARCH_CLASSID, "SNESLineSearch","Linesearch","SNESLineSearch",comm,SNESLineSearchDestroy,SNESLineSearchView);
194: linesearch->vec_sol_new = NULL;
195: linesearch->vec_func_new = NULL;
196: linesearch->vec_sol = NULL;
197: linesearch->vec_func = NULL;
198: linesearch->vec_update = NULL;
200: linesearch->lambda = 1.0;
201: linesearch->fnorm = 1.0;
202: linesearch->ynorm = 1.0;
203: linesearch->xnorm = 1.0;
204: linesearch->result = SNES_LINESEARCH_SUCCEEDED;
205: linesearch->norms = PETSC_TRUE;
206: linesearch->keeplambda = PETSC_FALSE;
207: linesearch->damping = 1.0;
208: linesearch->maxstep = 1e8;
209: linesearch->steptol = 1e-12;
210: linesearch->rtol = 1e-8;
211: linesearch->atol = 1e-15;
212: linesearch->ltol = 1e-8;
213: linesearch->precheckctx = NULL;
214: linesearch->postcheckctx = NULL;
215: linesearch->max_its = 1;
216: linesearch->setupcalled = PETSC_FALSE;
217: *outlinesearch = linesearch;
218: return(0);
219: }
221: /*@
222: SNESLineSearchSetUp - Prepares the line search for being applied by allocating
223: any required vectors.
225: Collective on SNESLineSearch
227: Input Parameters:
228: . linesearch - The LineSearch instance.
230: Notes:
231: For most cases, this needn't be called by users or outside of SNESLineSearchApply().
232: The only current case where this is called outside of this is for the VI
233: solvers, which modify the solution and work vectors before the first call
234: of SNESLineSearchApply, requiring the SNESLineSearch work vectors to be
235: allocated upfront.
237: Level: advanced
239: .keywords: SNESLineSearch, SetUp
241: .seealso: SNESGetLineSearch(), SNESLineSearchReset()
242: @*/
244: PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch)
245: {
249: if (!((PetscObject)linesearch)->type_name) {
250: SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);
251: }
252: if (!linesearch->setupcalled) {
253: if (!linesearch->vec_sol_new) {
254: VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new);
255: }
256: if (!linesearch->vec_func_new) {
257: VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new);
258: }
259: if (linesearch->ops->setup) {
260: (*linesearch->ops->setup)(linesearch);
261: }
262: if (!linesearch->ops->snesfunc) {SNESLineSearchSetFunction(linesearch,SNESComputeFunction);}
263: linesearch->lambda = linesearch->damping;
264: linesearch->setupcalled = PETSC_TRUE;
265: }
266: return(0);
267: }
270: /*@
271: SNESLineSearchReset - Undoes the SNESLineSearchSetUp() and deletes any Vecs or Mats allocated by the line search.
273: Collective on SNESLineSearch
275: Input Parameters:
276: . linesearch - The LineSearch instance.
278: Notes:
279: Usually only called by SNESReset()
281: Level: developer
283: .keywords: SNESLineSearch, Reset
285: .seealso: SNESGetLineSearch(), SNESLineSearchSetUp()
286: @*/
288: PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch)
289: {
293: if (linesearch->ops->reset) (*linesearch->ops->reset)(linesearch);
295: VecDestroy(&linesearch->vec_sol_new);
296: VecDestroy(&linesearch->vec_func_new);
298: VecDestroyVecs(linesearch->nwork, &linesearch->work);
300: linesearch->nwork = 0;
301: linesearch->setupcalled = PETSC_FALSE;
302: return(0);
303: }
305: /*@C
306: SNESLineSearchSetFunction - Sets the function evaluation used by the SNES line search
308: Input Parameters:
309: . linesearch - the SNESLineSearch context
310: + func - function evaluation routine
312: Level: developer
314: Notes:
315: This is used internally by PETSc and not called by users
317: .keywords: get, linesearch, pre-check
319: .seealso: SNESGetLineSearch(), SNESSetFunction()
320: @*/
321: PetscErrorCode SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES,Vec,Vec))
322: {
325: linesearch->ops->snesfunc = func;
326: return(0);
327: }
329: /*@C
330: SNESLineSearchSetPreCheck - Sets a user function that is called after the initial search direction has been computed but
331: before the line search routine has been applied. Allows the user to adjust the result of (usually a linear solve) that
332: determined the search direction.
334: Logically Collective on SNESLineSearch
336: Input Parameters:
337: + linesearch - the SNESLineSearch context
338: . func - [optional] function evaluation routine, see SNESLineSearchPreCheck() for the calling sequence
339: - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
341: Level: intermediate
343: .keywords: set, linesearch, pre-check
345: .seealso: SNESGetLineSearch(), SNESLineSearchPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
346: @*/
347: PetscErrorCode SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void *ctx)
348: {
351: if (func) linesearch->ops->precheck = func;
352: if (ctx) linesearch->precheckctx = ctx;
353: return(0);
354: }
356: /*@C
357: SNESLineSearchGetPreCheck - Gets the pre-check function for the line search routine.
359: Input Parameters:
360: . linesearch - the SNESLineSearch context
362: Output Parameters:
363: + func - [optional] function evaluation routine, see SNESLineSearchPreCheck() for calling sequence
364: - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
366: Level: intermediate
368: .keywords: get, linesearch, pre-check
370: .seealso: SNESGetLineSearch(), SNESGetLineSearch(), SNESLineSearchPreCheck(), SNESLineSearchGetPostCheck(), SNESLineSearchSetPreCheck(), SNESLineSearchSetPostCheck()
371: @*/
372: PetscErrorCode SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void **ctx)
373: {
376: if (func) *func = linesearch->ops->precheck;
377: if (ctx) *ctx = linesearch->precheckctx;
378: return(0);
379: }
381: /*@C
382: SNESLineSearchSetPostCheck - Sets a user function that is called after the line search has been applied to determine the step
383: direction and length. Allows the user a chance to change or override the decision of the line search routine
385: Logically Collective on SNESLineSearch
387: Input Parameters:
388: + linesearch - the SNESLineSearch context
389: . func - [optional] function evaluation routine, see SNESLineSearchPostCheck() for the calling sequence
390: - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
392: Level: intermediate
394: .keywords: set, linesearch, post-check
396: .seealso: SNESGetLineSearch(), SNESLineSearchPostCheck(), SNESLineSearchSetPreCheck(), SNESLineSearchGetPreCheck(), SNESLineSearchGetPostCheck()
397: @*/
398: PetscErrorCode SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void *ctx)
399: {
402: if (func) linesearch->ops->postcheck = func;
403: if (ctx) linesearch->postcheckctx = ctx;
404: return(0);
405: }
407: /*@C
408: SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine.
410: Input Parameters:
411: . linesearch - the SNESLineSearch context
413: Output Parameters:
414: + func - [optional] function evaluation routine, see for the calling sequence SNESLineSearchPostCheck()
415: - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
417: Level: intermediate
419: .keywords: get, linesearch, post-check
421: .seealso: SNESGetLineSearch(), SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchPostCheck(), SNESLineSearchSetPreCheck()
422: @*/
423: PetscErrorCode SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void **ctx)
424: {
427: if (func) *func = linesearch->ops->postcheck;
428: if (ctx) *ctx = linesearch->postcheckctx;
429: return(0);
430: }
432: /*@
433: SNESLineSearchPreCheck - Prepares the line search for being applied.
435: Logically Collective on SNESLineSearch
437: Input Parameters:
438: + linesearch - The linesearch instance.
439: . X - The current solution
440: - Y - The step direction
442: Output Parameters:
443: . changed - Indicator that the precheck routine has changed anything
445: Level: developer
447: .keywords: SNESLineSearch, Create
449: .seealso: SNESGetLineSearch(), SNESLineSearchPostCheck(), SNESLineSearchSetPreCheck(), SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck()
450: @*/
451: PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed)
452: {
456: *changed = PETSC_FALSE;
457: if (linesearch->ops->precheck) {
458: (*linesearch->ops->precheck)(linesearch, X, Y, changed, linesearch->precheckctx);
460: }
461: return(0);
462: }
464: /*@
465: SNESLineSearchPostCheck - Prepares the line search for being applied.
467: Logically Collective on SNESLineSearch
469: Input Parameters:
470: + linesearch - The linesearch context
471: . X - The last solution
472: . Y - The step direction
473: - W - The updated solution, W = X + lambda*Y for some lambda
475: Output Parameters:
476: + changed_Y - Indicator if the direction Y has been changed.
477: - changed_W - Indicator if the new candidate solution W has been changed.
479: Level: developer
481: .keywords: SNESLineSearch, Create
483: .seealso: SNESGetLineSearch(), SNESLineSearchPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck(), SNESLineSearchSetPrecheck(), SNESLineSearchGetPrecheck()
484: @*/
485: PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch,Vec X,Vec Y,Vec W,PetscBool *changed_Y,PetscBool *changed_W)
486: {
490: *changed_Y = PETSC_FALSE;
491: *changed_W = PETSC_FALSE;
492: if (linesearch->ops->postcheck) {
493: (*linesearch->ops->postcheck)(linesearch,X,Y,W,changed_Y,changed_W,linesearch->postcheckctx);
496: }
497: return(0);
498: }
500: /*@C
501: SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration
503: Logically Collective on SNESLineSearch
505: Input Arguments:
506: + linesearch - linesearch context
507: . X - base state for this step
508: . Y - initial correction
509: - ctx - context for this function
511: Output Arguments:
512: + Y - correction, possibly modified
513: - changed - flag indicating that Y was modified
515: Options Database Key:
516: + -snes_linesearch_precheck_picard - activate this routine
517: - -snes_linesearch_precheck_picard_angle - angle
519: Level: advanced
521: Notes:
522: This function should be passed to SNESLineSearchSetPreCheck()
524: The justification for this method involves the linear convergence of a Picard iteration
525: so the Picard linearization should be provided in place of the "Jacobian". This correction
526: is generally not useful when using a Newton linearization.
528: Reference:
529: Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology.
531: .seealso: SNESGetLineSearch(), SNESLineSearchSetPreCheck()
532: @*/
533: PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed,void *ctx)
534: {
536: PetscReal angle = *(PetscReal*)linesearch->precheckctx;
537: Vec Ylast;
538: PetscScalar dot;
539: PetscInt iter;
540: PetscReal ynorm,ylastnorm,theta,angle_radians;
541: SNES snes;
544: SNESLineSearchGetSNES(linesearch, &snes);
545: PetscObjectQuery((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject*)&Ylast);
546: if (!Ylast) {
547: VecDuplicate(Y,&Ylast);
548: PetscObjectCompose((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject)Ylast);
549: PetscObjectDereference((PetscObject)Ylast);
550: }
551: SNESGetIterationNumber(snes,&iter);
552: if (iter < 2) {
553: VecCopy(Y,Ylast);
554: *changed = PETSC_FALSE;
555: return(0);
556: }
558: VecDot(Y,Ylast,&dot);
559: VecNorm(Y,NORM_2,&ynorm);
560: VecNorm(Ylast,NORM_2,&ylastnorm);
561: /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
562: theta = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0));
563: angle_radians = angle * PETSC_PI / 180.;
564: if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
565: /* Modify the step Y */
566: PetscReal alpha,ydiffnorm;
567: VecAXPY(Ylast,-1.0,Y);
568: VecNorm(Ylast,NORM_2,&ydiffnorm);
569: alpha = (ydiffnorm > .001*ylastnorm) ? ylastnorm / ydiffnorm : 1000.0;
570: VecCopy(Y,Ylast);
571: VecScale(Y,alpha);
572: PetscInfo3(snes,"Angle %14.12e degrees less than threshold %14.12e, corrected step by alpha=%14.12e\n",(double)(theta*180./PETSC_PI),(double)angle,(double)alpha);
573: *changed = PETSC_TRUE;
574: } else {
575: PetscInfo2(snes,"Angle %14.12e degrees exceeds threshold %14.12e, no correction applied\n",(double)(theta*180./PETSC_PI),(double)angle);
576: VecCopy(Y,Ylast);
577: *changed = PETSC_FALSE;
578: }
579: return(0);
580: }
582: /*@
583: SNESLineSearchApply - Computes the line-search update.
585: Collective on SNESLineSearch
587: Input Parameters:
588: + linesearch - The linesearch context
589: . X - The current solution
590: . F - The current function
591: . fnorm - The current norm
592: - Y - The search direction
594: Output Parameters:
595: + X - The new solution
596: . F - The new function
597: - fnorm - The new function norm
599: Options Database Keys:
600: + -snes_linesearch_type - basic, bt, l2, cp, nleqerr, shell
601: . -snes_linesearch_monitor [:filename] - Print progress of line searches
602: . -snes_linesearch_damping - The linesearch damping parameter, default is 1.0 (no damping)
603: . -snes_linesearch_norms - Turn on/off the linesearch norms computation (SNESLineSearchSetComputeNorms())
604: . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess
605: - -snes_linesearch_max_it - The number of iterations for iterative line searches
607: Notes:
608: This is typically called from within a SNESSolve() implementation in order to
609: help with convergence of the nonlinear method. Various SNES types use line searches
610: in different ways, but the overarching theme is that a line search is used to determine
611: an optimal damping parameter of a step at each iteration of the method. Each
612: application of the line search may invoke SNESComputeFunction() several times, and
613: therefore may be fairly expensive.
615: Level: Intermediate
617: .keywords: SNESLineSearch, Create
619: .seealso: SNESGetLineSearch(), SNESLineSearchCreate(), SNESLineSearchPreCheck(), SNESLineSearchPostCheck(), SNESSolve(), SNESComputeFunction(), SNESLineSearchSetComputeNorms(),
620: SNESLineSearchType, SNESLineSearchSetType()
621: @*/
622: PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y)
623: {
632: linesearch->result = SNES_LINESEARCH_SUCCEEDED;
634: linesearch->vec_sol = X;
635: linesearch->vec_update = Y;
636: linesearch->vec_func = F;
638: SNESLineSearchSetUp(linesearch);
640: if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */
642: if (fnorm) linesearch->fnorm = *fnorm;
643: else {
644: VecNorm(F, NORM_2, &linesearch->fnorm);
645: }
647: PetscLogEventBegin(SNESLINESEARCH_Apply,linesearch,X,F,Y);
649: (*linesearch->ops->apply)(linesearch);
651: PetscLogEventEnd(SNESLINESEARCH_Apply,linesearch,X,F,Y);
653: if (fnorm) *fnorm = linesearch->fnorm;
654: return(0);
655: }
657: /*@
658: SNESLineSearchDestroy - Destroys the line search instance.
660: Collective on SNESLineSearch
662: Input Parameters:
663: . linesearch - The linesearch context
665: Level: developer
667: .keywords: SNESLineSearch, Destroy
669: .seealso: SNESGetLineSearch(), SNESLineSearchCreate(), SNESLineSearchReset(), SNESDestroy()
670: @*/
671: PetscErrorCode SNESLineSearchDestroy(SNESLineSearch * linesearch)
672: {
676: if (!*linesearch) return(0);
678: if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = 0; return(0);}
679: PetscObjectSAWsViewOff((PetscObject)*linesearch);
680: SNESLineSearchReset(*linesearch);
681: if ((*linesearch)->ops->destroy) (*linesearch)->ops->destroy(*linesearch);
682: PetscViewerDestroy(&(*linesearch)->monitor);
683: SNESLineSearchMonitorCancel((*linesearch));
684: PetscHeaderDestroy(linesearch);
685: return(0);
686: }
688: /*@
689: SNESLineSearchSetDefaultMonitor - Turns on/off printing useful information and debugging output about the line search.
691: Input Parameters:
692: + linesearch - the linesearch object
693: - viewer - an ASCII PetscViewer or NULL to turn off monitor
695: Logically Collective on SNESLineSearch
697: Options Database:
698: . -snes_linesearch_monitor [:filename] - enables the monitor
700: Level: intermediate
702: Developer Note: This monitor is implemented differently than the other SNESLineSearchMonitors that are set with
703: SNESLineSearchMonitorSet() since it is called in many locations of the line search routines to display aspects of the
704: line search that are not visible to the other monitors.
706: .seealso: SNESGetLineSearch(), SNESLineSearchGetDefaultMonitor(), PetscViewer, SNESLineSearchSetMonitor()
707: @*/
708: PetscErrorCode SNESLineSearchSetDefaultMonitor(SNESLineSearch linesearch, PetscViewer viewer)
709: {
713: if (viewer) {PetscObjectReference((PetscObject)viewer);}
714: PetscViewerDestroy(&linesearch->monitor);
715: linesearch->monitor = viewer;
716: return(0);
717: }
719: /*@
720: SNESLineSearchGetDefaultMonitor - Gets the PetscViewer instance for the line search monitor.
722: Input Parameter:
723: . linesearch - linesearch context
725: Output Parameter:
726: . monitor - monitor context
728: Logically Collective on SNES
730: Options Database Keys:
731: . -snes_linesearch_monitor - enables the monitor
733: Level: intermediate
735: .seealso: SNESGetLineSearch(), SNESLineSearchSetDefaultMonitor(), PetscViewer
736: @*/
737: PetscErrorCode SNESLineSearchGetDefaultMonitor(SNESLineSearch linesearch, PetscViewer *monitor)
738: {
741: if (monitor) {
743: *monitor = linesearch->monitor;
744: }
745: return(0);
746: }
748: /*@C
749: SNESLineSearchMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
751: Collective on SNESLineSearch
753: Input Parameters:
754: + ls - LineSearch object you wish to monitor
755: . name - the monitor type one is seeking
756: . help - message indicating what monitoring is done
757: . manual - manual page for the monitor
758: . monitor - the monitor function
759: - monitorsetup - a function that is called once ONLY if the user selected this monitor that may set additional features of the SNESLineSearch or PetscViewer objects
761: Level: developer
763: .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
764: PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
765: PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
766: PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
767: PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
768: PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
769: PetscOptionsFList(), PetscOptionsEList()
770: @*/
771: PetscErrorCode SNESLineSearchMonitorSetFromOptions(SNESLineSearch ls,const char name[],const char help[], const char manual[],PetscErrorCode (*monitor)(SNESLineSearch,PetscViewerAndFormat*),PetscErrorCode (*monitorsetup)(SNESLineSearch,PetscViewerAndFormat*))
772: {
773: PetscErrorCode ierr;
774: PetscViewer viewer;
775: PetscViewerFormat format;
776: PetscBool flg;
779: PetscOptionsGetViewer(PetscObjectComm((PetscObject)ls),((PetscObject) ls)->options,((PetscObject)ls)->prefix,name,&viewer,&format,&flg);
780: if (flg) {
781: PetscViewerAndFormat *vf;
782: PetscViewerAndFormatCreate(viewer,format,&vf);
783: PetscObjectDereference((PetscObject)viewer);
784: if (monitorsetup) {
785: (*monitorsetup)(ls,vf);
786: }
787: SNESLineSearchMonitorSet(ls,(PetscErrorCode (*)(SNESLineSearch,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);
788: }
789: return(0);
790: }
792: /*@
793: SNESLineSearchSetFromOptions - Sets options for the line search
795: Input Parameters:
796: . linesearch - linesearch context
798: Options Database Keys:
799: + -snes_linesearch_type <type> - basic, bt, l2, cp, nleqerr, shell
800: . -snes_linesearch_order <order> - 1, 2, 3. Most types only support certain orders (bt supports 2 or 3)
801: . -snes_linesearch_norms - Turn on/off the linesearch norms for the basic linesearch typem (SNESLineSearchSetComputeNorms())
802: . -snes_linesearch_minlambda - The minimum step length
803: . -snes_linesearch_maxstep - The maximum step size
804: . -snes_linesearch_rtol - Relative tolerance for iterative line searches
805: . -snes_linesearch_atol - Absolute tolerance for iterative line searches
806: . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches
807: . -snes_linesearch_max_it - The number of iterations for iterative line searches
808: . -snes_linesearch_monitor [:filename] - Print progress of line searches
809: . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine
810: . -snes_linesearch_damping - The linesearch damping parameter
811: . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess.
812: . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method
813: - -snes_linesearch_precheck_picard_angle - Angle used in picard precheck method
815: Logically Collective on SNESLineSearch
817: Level: intermediate
819: .seealso: SNESLineSearchCreate(), SNESLineSearchSetOrder(), SNESLineSearchSetType(), SNESLineSearchSetTolerances(), SNESLineSearchSetDamping(), SNESLineSearchPreCheckPicard(),
820: SNESLineSearchType, SNESLineSearchSetComputeNorms()
821: @*/
822: PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch)
823: {
824: PetscErrorCode ierr;
825: const char *deft = SNESLINESEARCHBASIC;
826: char type[256];
827: PetscBool flg, set;
828: PetscViewer viewer;
831: SNESLineSearchRegisterAll();
833: PetscObjectOptionsBegin((PetscObject)linesearch);
834: if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name;
835: PetscOptionsFList("-snes_linesearch_type","Linesearch type","SNESLineSearchSetType",SNESLineSearchList,deft,type,256,&flg);
836: if (flg) {
837: SNESLineSearchSetType(linesearch,type);
838: } else if (!((PetscObject)linesearch)->type_name) {
839: SNESLineSearchSetType(linesearch,deft);
840: }
842: PetscOptionsGetViewer(PetscObjectComm((PetscObject)linesearch),((PetscObject) linesearch)->options,((PetscObject)linesearch)->prefix,"-snes_linesearch_monitor",&viewer,NULL,&set);
843: if (set) {
844: SNESLineSearchSetDefaultMonitor(linesearch,viewer);
845: PetscViewerDestroy(&viewer);
846: }
847: SNESLineSearchMonitorSetFromOptions(linesearch,"-snes_linesearch_monitor_solution_update","View correction at each iteration","SNESLineSearchMonitorSolutionUpdate",SNESLineSearchMonitorSolutionUpdate,NULL);
849: /* tolerances */
850: PetscOptionsReal("-snes_linesearch_minlambda","Minimum step length","SNESLineSearchSetTolerances",linesearch->steptol,&linesearch->steptol,NULL);
851: PetscOptionsReal("-snes_linesearch_maxstep","Maximum step size","SNESLineSearchSetTolerances",linesearch->maxstep,&linesearch->maxstep,NULL);
852: PetscOptionsReal("-snes_linesearch_rtol","Relative tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->rtol,&linesearch->rtol,NULL);
853: PetscOptionsReal("-snes_linesearch_atol","Absolute tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->atol,&linesearch->atol,NULL);
854: PetscOptionsReal("-snes_linesearch_ltol","Change in lambda tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->ltol,&linesearch->ltol,NULL);
855: PetscOptionsInt("-snes_linesearch_max_it","Maximum iterations for iterative line searches","SNESLineSearchSetTolerances",linesearch->max_its,&linesearch->max_its,NULL);
857: /* damping parameters */
858: PetscOptionsReal("-snes_linesearch_damping","Line search damping and initial step guess","SNESLineSearchSetDamping",linesearch->damping,&linesearch->damping,NULL);
860: PetscOptionsBool("-snes_linesearch_keeplambda","Use previous lambda as damping","SNESLineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,NULL);
862: /* precheck */
863: PetscOptionsBool("-snes_linesearch_precheck_picard","Use a correction that sometimes improves convergence of Picard iteration","SNESLineSearchPreCheckPicard",flg,&flg,&set);
864: if (set) {
865: if (flg) {
866: linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */
868: PetscOptionsReal("-snes_linesearch_precheck_picard_angle","Maximum angle at which to activate the correction",
869: "none",linesearch->precheck_picard_angle,&linesearch->precheck_picard_angle,NULL);
870: SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&linesearch->precheck_picard_angle);
871: } else {
872: SNESLineSearchSetPreCheck(linesearch,NULL,NULL);
873: }
874: }
875: PetscOptionsInt("-snes_linesearch_order","Order of approximation used in the line search","SNESLineSearchSetOrder",linesearch->order,&linesearch->order,NULL);
876: PetscOptionsBool("-snes_linesearch_norms","Compute final norms in line search","SNESLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,NULL);
878: if (linesearch->ops->setfromoptions) {
879: (*linesearch->ops->setfromoptions)(PetscOptionsObject,linesearch);
880: }
882: PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)linesearch);
883: PetscOptionsEnd();
884: return(0);
885: }
887: /*@
888: SNESLineSearchView - Prints useful information about the line search
890: Input Parameters:
891: . linesearch - linesearch context
893: Logically Collective on SNESLineSearch
895: Level: intermediate
897: .seealso: SNESLineSearchCreate()
898: @*/
899: PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer)
900: {
902: PetscBool iascii;
906: if (!viewer) {
907: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch),&viewer);
908: }
912: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
913: if (iascii) {
914: PetscObjectPrintClassNamePrefixType((PetscObject)linesearch,viewer);
915: if (linesearch->ops->view) {
916: PetscViewerASCIIPushTab(viewer);
917: (*linesearch->ops->view)(linesearch,viewer);
918: PetscViewerASCIIPopTab(viewer);
919: }
920: PetscViewerASCIIPrintf(viewer," maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep,(double)linesearch->steptol);
921: PetscViewerASCIIPrintf(viewer," tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol,(double)linesearch->atol,(double)linesearch->ltol);
922: PetscViewerASCIIPrintf(viewer," maximum iterations=%D\n", linesearch->max_its);
923: if (linesearch->ops->precheck) {
924: if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) {
925: PetscViewerASCIIPrintf(viewer," using precheck step to speed up Picard convergence\n", linesearch->max_its);
926: } else {
927: PetscViewerASCIIPrintf(viewer," using user-defined precheck step\n", linesearch->max_its);
928: }
929: }
930: if (linesearch->ops->postcheck) {
931: PetscViewerASCIIPrintf(viewer," using user-defined postcheck step\n", linesearch->max_its);
932: }
933: }
934: return(0);
935: }
937: /*@C
938: SNESLineSearchSetType - Sets the linesearch type
940: Logically Collective on SNESLineSearch
942: Input Parameters:
943: + linesearch - linesearch context
944: - type - The type of line search to be used
946: Available Types:
947: + SNESLINESEARCHBASIC - Simple damping line search, defaults to using the full Newton step
948: . SNESLINESEARCHBT - Backtracking line search over the L2 norm of the function
949: . SNESLINESEARCHL2 - Secant line search over the L2 norm of the function
950: . SNESLINESEARCHCP - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x)
951: . SNESLINESEARCHNLEQERR - Affine-covariant error-oriented linesearch
952: - SNESLINESEARCHSHELL - User provided SNESLineSearch implementation
954: Options Database:
955: . -snes_linesearch_type <type> - basic, bt, l2, cp, nleqerr, shell
957: Level: intermediate
959: .seealso: SNESLineSearchCreate(), SNESLineSearchType, SNESLineSearchSetFromOptions()
960: @*/
961: PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type)
962: {
963: PetscErrorCode ierr,(*r)(SNESLineSearch);
964: PetscBool match;
970: PetscObjectTypeCompare((PetscObject)linesearch,type,&match);
971: if (match) return(0);
973: PetscFunctionListFind(SNESLineSearchList,type,&r);
974: if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type);
975: /* Destroy the previous private linesearch context */
976: if (linesearch->ops->destroy) {
977: (*(linesearch)->ops->destroy)(linesearch);
979: linesearch->ops->destroy = NULL;
980: }
981: /* Reinitialize function pointers in SNESLineSearchOps structure */
982: linesearch->ops->apply = 0;
983: linesearch->ops->view = 0;
984: linesearch->ops->setfromoptions = 0;
985: linesearch->ops->destroy = 0;
987: PetscObjectChangeTypeName((PetscObject)linesearch,type);
988: (*r)(linesearch);
989: return(0);
990: }
992: /*@
993: SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation.
995: Input Parameters:
996: + linesearch - linesearch context
997: - snes - The snes instance
999: Level: developer
1001: Notes:
1002: This happens automatically when the line search is obtained/created with
1003: SNESGetLineSearch(). This routine is therefore mainly called within SNES
1004: implementations.
1006: Level: developer
1008: .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
1009: @*/
1010: PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes)
1011: {
1015: linesearch->snes = snes;
1016: return(0);
1017: }
1019: /*@
1020: SNESLineSearchGetSNES - Gets the SNES instance associated with the line search.
1021: Having an associated SNES is necessary because most line search implementations must be able to
1022: evaluate the function using SNESComputeFunction() for the associated SNES. This routine
1023: is used in the line search implementations when one must get this associated SNES instance.
1025: Input Parameters:
1026: . linesearch - linesearch context
1028: Output Parameters:
1029: . snes - The snes instance
1031: Level: developer
1033: .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
1034: @*/
1035: PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes)
1036: {
1040: *snes = linesearch->snes;
1041: return(0);
1042: }
1044: /*@
1045: SNESLineSearchGetLambda - Gets the last linesearch steplength discovered.
1047: Input Parameters:
1048: . linesearch - linesearch context
1050: Output Parameters:
1051: . lambda - The last steplength computed during SNESLineSearchApply()
1053: Level: advanced
1055: Notes:
1056: This is useful in methods where the solver is ill-scaled and
1057: requires some adaptive notion of the difference in scale between the
1058: solution and the function. For instance, SNESQN may be scaled by the
1059: line search lambda using the argument -snes_qn_scaling ls.
1061: .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetDamping(), SNESLineSearchApply()
1062: @*/
1063: PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda)
1064: {
1068: *lambda = linesearch->lambda;
1069: return(0);
1070: }
1072: /*@
1073: SNESLineSearchSetLambda - Sets the linesearch steplength.
1075: Input Parameters:
1076: + linesearch - linesearch context
1077: - lambda - The last steplength.
1079: Notes:
1080: This routine is typically used within implementations of SNESLineSearchApply()
1081: to set the final steplength. This routine (and SNESLineSearchGetLambda()) were
1082: added in order to facilitate Quasi-Newton methods that use the previous steplength
1083: as an inner scaling parameter.
1085: Level: advanced
1087: .seealso: SNESLineSearchGetLambda()
1088: @*/
1089: PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
1090: {
1093: linesearch->lambda = lambda;
1094: return(0);
1095: }
1097: /*@
1098: SNESLineSearchGetTolerances - Gets the tolerances for the linesearch. These include
1099: tolerances for the relative and absolute change in the function norm, the change
1100: in lambda for iterative line searches, the minimum steplength, the maximum steplength,
1101: and the maximum number of iterations the line search procedure may take.
1103: Input Parameters:
1104: . linesearch - linesearch context
1106: Output Parameters:
1107: + steptol - The minimum steplength
1108: . maxstep - The maximum steplength
1109: . rtol - The relative tolerance for iterative line searches
1110: . atol - The absolute tolerance for iterative line searches
1111: . ltol - The change in lambda tolerance for iterative line searches
1112: - max_it - The maximum number of iterations of the line search
1114: Level: intermediate
1116: Notes:
1117: Different line searches may implement these parameters slightly differently as
1118: the type requires.
1120: .seealso: SNESLineSearchSetTolerances()
1121: @*/
1122: PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
1123: {
1126: if (steptol) {
1128: *steptol = linesearch->steptol;
1129: }
1130: if (maxstep) {
1132: *maxstep = linesearch->maxstep;
1133: }
1134: if (rtol) {
1136: *rtol = linesearch->rtol;
1137: }
1138: if (atol) {
1140: *atol = linesearch->atol;
1141: }
1142: if (ltol) {
1144: *ltol = linesearch->ltol;
1145: }
1146: if (max_its) {
1148: *max_its = linesearch->max_its;
1149: }
1150: return(0);
1151: }
1153: /*@
1154: SNESLineSearchSetTolerances - Gets the tolerances for the linesearch. These include
1155: tolerances for the relative and absolute change in the function norm, the change
1156: in lambda for iterative line searches, the minimum steplength, the maximum steplength,
1157: and the maximum number of iterations the line search procedure may take.
1159: Input Parameters:
1160: + linesearch - linesearch context
1161: . steptol - The minimum steplength
1162: . maxstep - The maximum steplength
1163: . rtol - The relative tolerance for iterative line searches
1164: . atol - The absolute tolerance for iterative line searches
1165: . ltol - The change in lambda tolerance for iterative line searches
1166: - max_it - The maximum number of iterations of the line search
1168: Notes:
1169: The user may choose to not set any of the tolerances using PETSC_DEFAULT in
1170: place of an argument.
1172: Level: intermediate
1174: .seealso: SNESLineSearchGetTolerances()
1175: @*/
1176: PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its)
1177: {
1187: if (steptol!= PETSC_DEFAULT) {
1188: if (steptol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Minimum step length %14.12e must be non-negative",(double)steptol);
1189: linesearch->steptol = steptol;
1190: }
1192: if (maxstep!= PETSC_DEFAULT) {
1193: if (maxstep < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum step length %14.12e must be non-negative",(double)maxstep);
1194: linesearch->maxstep = maxstep;
1195: }
1197: if (rtol != PETSC_DEFAULT) {
1198: if (rtol < 0.0 || 1.0 <= rtol) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Relative tolerance %14.12e must be non-negative and less than 1.0",(double)rtol);
1199: linesearch->rtol = rtol;
1200: }
1202: if (atol != PETSC_DEFAULT) {
1203: if (atol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %14.12e must be non-negative",(double)atol);
1204: linesearch->atol = atol;
1205: }
1207: if (ltol != PETSC_DEFAULT) {
1208: if (ltol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Labmda tolerance %14.12e must be non-negative",(double)ltol);
1209: linesearch->ltol = ltol;
1210: }
1212: if (max_its != PETSC_DEFAULT) {
1213: if (max_its < 0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",max_its);
1214: linesearch->max_its = max_its;
1215: }
1216: return(0);
1217: }
1219: /*@
1220: SNESLineSearchGetDamping - Gets the line search damping parameter.
1222: Input Parameters:
1223: . linesearch - linesearch context
1225: Output Parameters:
1226: . damping - The damping parameter
1228: Level: advanced
1230: .seealso: SNESLineSearchGetStepTolerance(), SNESQN
1231: @*/
1233: PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping)
1234: {
1238: *damping = linesearch->damping;
1239: return(0);
1240: }
1242: /*@
1243: SNESLineSearchSetDamping - Sets the line search damping paramter.
1245: Input Parameters:
1246: + linesearch - linesearch context
1247: - damping - The damping parameter
1249: Options Database:
1250: . -snes_linesearch_damping
1251: Level: intermediate
1253: Notes:
1254: The basic line search merely takes the update step scaled by the damping parameter.
1255: The use of the damping parameter in the l2 and cp line searches is much more subtle;
1256: it is used as a starting point in calculating the secant step. However, the eventual
1257: step may be of greater length than the damping parameter. In the bt line search it is
1258: used as the maximum possible step length, as the bt line search only backtracks.
1260: .seealso: SNESLineSearchGetDamping()
1261: @*/
1262: PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping)
1263: {
1266: linesearch->damping = damping;
1267: return(0);
1268: }
1270: /*@
1271: SNESLineSearchGetOrder - Gets the line search approximation order.
1273: Input Parameters:
1274: . linesearch - linesearch context
1276: Output Parameters:
1277: . order - The order
1279: Possible Values for order:
1280: + 1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
1281: . 2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
1282: - 3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
1284: Level: intermediate
1286: .seealso: SNESLineSearchSetOrder()
1287: @*/
1289: PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch linesearch,PetscInt *order)
1290: {
1294: *order = linesearch->order;
1295: return(0);
1296: }
1298: /*@
1299: SNESLineSearchSetOrder - Sets the maximum order of the polynomial fit used in the line search
1301: Input Parameters:
1302: . linesearch - linesearch context
1303: . order - The damping parameter
1305: Level: intermediate
1307: Possible Values for order:
1308: + 1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
1309: . 2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
1310: - 3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
1312: Notes:
1313: Variable orders are supported by the following line searches:
1314: + bt - cubic and quadratic
1315: - cp - linear and quadratic
1317: .seealso: SNESLineSearchGetOrder(), SNESLineSearchSetDamping()
1318: @*/
1319: PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch linesearch,PetscInt order)
1320: {
1323: linesearch->order = order;
1324: return(0);
1325: }
1327: /*@
1328: SNESLineSearchGetNorms - Gets the norms for for X, Y, and F.
1330: Input Parameters:
1331: . linesearch - linesearch context
1333: Output Parameters:
1334: + xnorm - The norm of the current solution
1335: . fnorm - The norm of the current function
1336: - ynorm - The norm of the current update
1338: Notes:
1339: This function is mainly called from SNES implementations.
1341: Level: developer
1343: .seealso: SNESLineSearchSetNorms() SNESLineSearchGetVecs()
1344: @*/
1345: PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm)
1346: {
1349: if (xnorm) *xnorm = linesearch->xnorm;
1350: if (fnorm) *fnorm = linesearch->fnorm;
1351: if (ynorm) *ynorm = linesearch->ynorm;
1352: return(0);
1353: }
1355: /*@
1356: SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F.
1358: Input Parameters:
1359: + linesearch - linesearch context
1360: . xnorm - The norm of the current solution
1361: . fnorm - The norm of the current function
1362: - ynorm - The norm of the current update
1364: Level: advanced
1366: .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1367: @*/
1368: PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
1369: {
1372: linesearch->xnorm = xnorm;
1373: linesearch->fnorm = fnorm;
1374: linesearch->ynorm = ynorm;
1375: return(0);
1376: }
1378: /*@
1379: SNESLineSearchComputeNorms - Computes the norms of X, F, and Y.
1381: Input Parameters:
1382: . linesearch - linesearch context
1384: Options Database Keys:
1385: . -snes_linesearch_norms - turn norm computation on or off
1387: Level: intermediate
1389: .seealso: SNESLineSearchGetNorms, SNESLineSearchSetNorms(), SNESLineSearchSetComputeNorms()
1390: @*/
1391: PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
1392: {
1394: SNES snes;
1397: if (linesearch->norms) {
1398: if (linesearch->ops->vinorm) {
1399: SNESLineSearchGetSNES(linesearch, &snes);
1400: VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);
1401: VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);
1402: (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);
1403: } else {
1404: VecNormBegin(linesearch->vec_func, NORM_2, &linesearch->fnorm);
1405: VecNormBegin(linesearch->vec_sol, NORM_2, &linesearch->xnorm);
1406: VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);
1407: VecNormEnd(linesearch->vec_func, NORM_2, &linesearch->fnorm);
1408: VecNormEnd(linesearch->vec_sol, NORM_2, &linesearch->xnorm);
1409: VecNormEnd(linesearch->vec_update, NORM_2, &linesearch->ynorm);
1410: }
1411: }
1412: return(0);
1413: }
1415: /*@
1416: SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search.
1418: Input Parameters:
1419: + linesearch - linesearch context
1420: - flg - indicates whether or not to compute norms
1422: Options Database Keys:
1423: . -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic linesearch
1425: Notes:
1426: This is most relevant to the SNESLINESEARCHBASIC line search type since most line searches have a stopping criteria involving the norm.
1428: Level: intermediate
1430: .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetNorms(), SNESLineSearchComputeNorms(), SNESLINESEARCHBASIC
1431: @*/
1432: PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg)
1433: {
1435: linesearch->norms = flg;
1436: return(0);
1437: }
1439: /*@
1440: SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context
1442: Input Parameters:
1443: . linesearch - linesearch context
1445: Output Parameters:
1446: + X - Solution vector
1447: . F - Function vector
1448: . Y - Search direction vector
1449: . W - Solution work vector
1450: - G - Function work vector
1452: Notes:
1453: At the beginning of a line search application, X should contain a
1454: solution and the vector F the function computed at X. At the end of the
1455: line search application, X should contain the new solution, and F the
1456: function evaluated at the new solution.
1458: These vectors are owned by the SNESLineSearch and should not be destroyed by the caller
1460: Level: advanced
1462: .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1463: @*/
1464: PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G)
1465: {
1468: if (X) {
1470: *X = linesearch->vec_sol;
1471: }
1472: if (F) {
1474: *F = linesearch->vec_func;
1475: }
1476: if (Y) {
1478: *Y = linesearch->vec_update;
1479: }
1480: if (W) {
1482: *W = linesearch->vec_sol_new;
1483: }
1484: if (G) {
1486: *G = linesearch->vec_func_new;
1487: }
1488: return(0);
1489: }
1491: /*@
1492: SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context
1494: Input Parameters:
1495: + linesearch - linesearch context
1496: . X - Solution vector
1497: . F - Function vector
1498: . Y - Search direction vector
1499: . W - Solution work vector
1500: - G - Function work vector
1502: Level: advanced
1504: .seealso: SNESLineSearchSetNorms(), SNESLineSearchGetVecs()
1505: @*/
1506: PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G)
1507: {
1510: if (X) {
1512: linesearch->vec_sol = X;
1513: }
1514: if (F) {
1516: linesearch->vec_func = F;
1517: }
1518: if (Y) {
1520: linesearch->vec_update = Y;
1521: }
1522: if (W) {
1524: linesearch->vec_sol_new = W;
1525: }
1526: if (G) {
1528: linesearch->vec_func_new = G;
1529: }
1530: return(0);
1531: }
1533: /*@C
1534: SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1535: SNES options in the database.
1537: Logically Collective on SNESLineSearch
1539: Input Parameters:
1540: + snes - the SNES context
1541: - prefix - the prefix to prepend to all option names
1543: Notes:
1544: A hyphen (-) must NOT be given at the beginning of the prefix name.
1545: The first character of all runtime options is AUTOMATICALLY the hyphen.
1547: Level: advanced
1549: .keywords: SNESLineSearch, append, options, prefix, database
1551: .seealso: SNESGetOptionsPrefix()
1552: @*/
1553: PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[])
1554: {
1559: PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);
1560: return(0);
1561: }
1563: /*@C
1564: SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all
1565: SNESLineSearch options in the database.
1567: Not Collective
1569: Input Parameter:
1570: . linesearch - the SNESLineSearch context
1572: Output Parameter:
1573: . prefix - pointer to the prefix string used
1575: Notes:
1576: On the fortran side, the user should pass in a string 'prefix' of
1577: sufficient length to hold the prefix.
1579: Level: advanced
1581: .keywords: SNESLineSearch, get, options, prefix, database
1583: .seealso: SNESAppendOptionsPrefix()
1584: @*/
1585: PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[])
1586: {
1591: PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);
1592: return(0);
1593: }
1595: /*@C
1596: SNESLineSearchSetWorkVecs - Gets work vectors for the line search.
1598: Input Parameter:
1599: + linesearch - the SNESLineSearch context
1600: - nwork - the number of work vectors
1602: Level: developer
1604: .keywords: SNESLineSearch, work, vector
1606: .seealso: SNESSetWorkVecs()
1607: @*/
1608: PetscErrorCode SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork)
1609: {
1613: if (linesearch->vec_sol) {
1614: VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);
1615: } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
1616: return(0);
1617: }
1619: /*@
1620: SNESLineSearchGetReason - Gets the success/failure status of the last line search application
1622: Input Parameters:
1623: . linesearch - linesearch context
1625: Output Parameters:
1626: . result - The success or failure status
1628: Notes:
1629: This is typically called after SNESLineSearchApply() in order to determine if the line-search failed
1630: (and set the SNES convergence accordingly).
1632: Level: intermediate
1634: .seealso: SNESLineSearchSetReason(), SNESLineSearchReason
1635: @*/
1636: PetscErrorCode SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result)
1637: {
1641: *result = linesearch->result;
1642: return(0);
1643: }
1645: /*@
1646: SNESLineSearchSetReason - Sets the success/failure status of the last line search application
1648: Input Parameters:
1649: + linesearch - linesearch context
1650: - result - The success or failure status
1652: Notes:
1653: This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set
1654: the success or failure of the line search method.
1656: Level: developer
1658: .seealso: SNESLineSearchGetSResult()
1659: @*/
1660: PetscErrorCode SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result)
1661: {
1664: linesearch->result = result;
1665: return(0);
1666: }
1668: /*@C
1669: SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.
1671: Input Parameters:
1672: + snes - nonlinear context obtained from SNESCreate()
1673: . projectfunc - function for projecting the function to the bounds
1674: - normfunc - function for computing the norm of an active set
1676: Logically Collective on SNES
1678: Calling sequence of projectfunc:
1679: .vb
1680: projectfunc (SNES snes, Vec X)
1681: .ve
1683: Input parameters for projectfunc:
1684: + snes - nonlinear context
1685: - X - current solution
1687: Output parameters for projectfunc:
1688: . X - Projected solution
1690: Calling sequence of normfunc:
1691: .vb
1692: projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm)
1693: .ve
1695: Input parameters for normfunc:
1696: + snes - nonlinear context
1697: . X - current solution
1698: - F - current residual
1700: Output parameters for normfunc:
1701: . fnorm - VI-specific norm of the function
1703: Notes:
1704: The VI solvers require projection of the solution to the feasible set. projectfunc should implement this.
1706: The VI solvers require special evaluation of the function norm such that the norm is only calculated
1707: on the inactive set. This should be implemented by normfunc.
1709: Level: developer
1711: .keywords: SNES, line search, VI, nonlinear, set, line search
1713: .seealso: SNESLineSearchGetVIFunctions(), SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck()
1714: @*/
1715: PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc)
1716: {
1719: if (projectfunc) linesearch->ops->viproject = projectfunc;
1720: if (normfunc) linesearch->ops->vinorm = normfunc;
1721: return(0);
1722: }
1724: /*@C
1725: SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.
1727: Input Parameters:
1728: . linesearch - the line search context, obtain with SNESGetLineSearch()
1730: Output Parameters:
1731: + projectfunc - function for projecting the function to the bounds
1732: - normfunc - function for computing the norm of an active set
1734: Logically Collective on SNES
1736: Level: developer
1738: .keywords: SNES, line search, VI, nonlinear, get, line search
1740: .seealso: SNESLineSearchSetVIFunctions(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
1741: @*/
1742: PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc)
1743: {
1745: if (projectfunc) *projectfunc = linesearch->ops->viproject;
1746: if (normfunc) *normfunc = linesearch->ops->vinorm;
1747: return(0);
1748: }
1750: /*@C
1751: SNESLineSearchRegister - See SNESLineSearchRegister()
1753: Level: advanced
1754: @*/
1755: PetscErrorCode SNESLineSearchRegister(const char sname[],PetscErrorCode (*function)(SNESLineSearch))
1756: {
1760: SNESInitializePackage();
1761: PetscFunctionListAdd(&SNESLineSearchList,sname,function);
1762: return(0);
1763: }