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 it
 28:   that one.

 30: .seealso: [](ch_snes), `SNES`, `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 `SNESLineSearch` implementations.
 57:   It does not typically need to be called by the user.

 59: .seealso: [](ch_snes), `SNES`, `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 monitor routine (use `NULL` if no context is desired)
 81: - monitordestroy - [optional] routine that frees monitor context (may be `NULL`)

 83:   Calling sequence of `f`:
 84: + ls   - the `SNESLineSearch` context
 85: - mctx - [optional] user-defined context for private data for the monitor routine

 87:   Calling sequence of `monitordestroy`:
 88: . mctx - [optional] user-defined context for private data for the monitor routine

 90:   Level: intermediate

 92:   Note:
 93:   Several different monitoring routines may be set by calling
 94:   `SNESLineSearchMonitorSet()` multiple times; all will be called in the
 95:   order in which they were set.

 97:   Fortran Note:
 98:   Only a single monitor function can be set for each `SNESLineSearch` object

100: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorDefault()`, `SNESLineSearchMonitorCancel()`
101: @*/
102: PetscErrorCode SNESLineSearchMonitorSet(SNESLineSearch ls, PetscErrorCode (*f)(SNESLineSearch ls, void *mctx), void *mctx, PetscErrorCode (*monitordestroy)(void **mctx))
103: {
104:   PetscInt  i;
105:   PetscBool identical;

107:   PetscFunctionBegin;
109:   for (i = 0; i < ls->numbermonitors; i++) {
110:     PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))f, mctx, monitordestroy, (PetscErrorCode(*)(void))ls->monitorftns[i], ls->monitorcontext[i], ls->monitordestroy[i], &identical));
111:     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
112:   }
113:   PetscCheck(ls->numbermonitors < MAXSNESLSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many monitors set");
114:   ls->monitorftns[ls->numbermonitors]      = f;
115:   ls->monitordestroy[ls->numbermonitors]   = monitordestroy;
116:   ls->monitorcontext[ls->numbermonitors++] = (void *)mctx;
117:   PetscFunctionReturn(PETSC_SUCCESS);
118: }

120: /*@C
121:   SNESLineSearchMonitorSolutionUpdate - Monitors each update of the function value the linesearch tries

123:   Collective

125:   Input Parameters:
126: + ls - the `SNESLineSearch` object
127: - vf - the context for the monitor, in this case it is an `PetscViewerAndFormat`

129:   Options Database Key:
130: . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine

132:   Level: developer

134:   This is not normally called directly but is passed to `SNESLineSearchMonitorSet()`

136: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorSet()`, `SNESMonitorSolution()`
137: @*/
138: PetscErrorCode SNESLineSearchMonitorSolutionUpdate(SNESLineSearch ls, PetscViewerAndFormat *vf)
139: {
140:   PetscViewer viewer = vf->viewer;
141:   Vec         Y, W, G;

143:   PetscFunctionBegin;
144:   PetscCall(SNESLineSearchGetVecs(ls, NULL, NULL, &Y, &W, &G));
145:   PetscCall(PetscViewerPushFormat(viewer, vf->format));
146:   PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted update to solution \n"));
147:   PetscCall(VecView(Y, viewer));
148:   PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted new solution \n"));
149:   PetscCall(VecView(W, viewer));
150:   PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted updated function value\n"));
151:   PetscCall(VecView(G, viewer));
152:   PetscCall(PetscViewerPopFormat(viewer));
153:   PetscFunctionReturn(PETSC_SUCCESS);
154: }

156: /*@
157:   SNESLineSearchCreate - Creates a `SNESLineSearch` context.

159:   Logically Collective

161:   Input Parameter:
162: . comm - MPI communicator for the line search (typically from the associated `SNES` context).

164:   Output Parameter:
165: . outlinesearch - the new line search context

167:   Level: developer

169:   Note:
170:   The preferred calling sequence is to use `SNESGetLineSearch()` to acquire the `SNESLineSearch` instance
171:   already associated with the `SNES`.

173: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `LineSearchDestroy()`, `SNESGetLineSearch()`
174: @*/
175: PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch)
176: {
177:   SNESLineSearch linesearch;

179:   PetscFunctionBegin;
180:   PetscAssertPointer(outlinesearch, 2);
181:   PetscCall(SNESInitializePackage());
182:   *outlinesearch = NULL;

184:   PetscCall(PetscHeaderCreate(linesearch, SNESLINESEARCH_CLASSID, "SNESLineSearch", "Linesearch", "SNESLineSearch", comm, SNESLineSearchDestroy, SNESLineSearchView));

186:   linesearch->vec_sol_new  = NULL;
187:   linesearch->vec_func_new = NULL;
188:   linesearch->vec_sol      = NULL;
189:   linesearch->vec_func     = NULL;
190:   linesearch->vec_update   = NULL;

192:   linesearch->lambda       = 1.0;
193:   linesearch->fnorm        = 1.0;
194:   linesearch->ynorm        = 1.0;
195:   linesearch->xnorm        = 1.0;
196:   linesearch->result       = SNES_LINESEARCH_SUCCEEDED;
197:   linesearch->norms        = PETSC_TRUE;
198:   linesearch->keeplambda   = PETSC_FALSE;
199:   linesearch->damping      = 1.0;
200:   linesearch->maxstep      = 1e8;
201:   linesearch->steptol      = 1e-12;
202:   linesearch->rtol         = 1e-8;
203:   linesearch->atol         = 1e-15;
204:   linesearch->ltol         = 1e-8;
205:   linesearch->precheckctx  = NULL;
206:   linesearch->postcheckctx = NULL;
207:   linesearch->max_its      = 1;
208:   linesearch->setupcalled  = PETSC_FALSE;
209:   linesearch->monitor      = NULL;
210:   *outlinesearch           = linesearch;
211:   PetscFunctionReturn(PETSC_SUCCESS);
212: }

214: /*@
215:   SNESLineSearchSetUp - Prepares the line search for being applied by allocating
216:   any required vectors.

218:   Collective

220:   Input Parameter:
221: . linesearch - The `SNESLineSearch` instance.

223:   Level: advanced

225:   Note:
226:   For most cases, this needn't be called by users or outside of `SNESLineSearchApply()`.
227:   The only current case where this is called outside of this is for the VI
228:   solvers, which modify the solution and work vectors before the first call
229:   of `SNESLineSearchApply()`, requiring the `SNESLineSearch` work vectors to be
230:   allocated upfront.

232: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchReset()`
233: @*/
234: PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch)
235: {
236:   PetscFunctionBegin;
237:   if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC));
238:   if (!linesearch->setupcalled) {
239:     if (!linesearch->vec_sol_new) PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new));
240:     if (!linesearch->vec_func_new) PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new));
241:     if (linesearch->ops->setup) PetscUseTypeMethod(linesearch, setup);
242:     if (!linesearch->ops->snesfunc) PetscCall(SNESLineSearchSetFunction(linesearch, SNESComputeFunction));
243:     linesearch->lambda      = linesearch->damping;
244:     linesearch->setupcalled = PETSC_TRUE;
245:   }
246:   PetscFunctionReturn(PETSC_SUCCESS);
247: }

249: /*@
250:   SNESLineSearchReset - Undoes the `SNESLineSearchSetUp()` and deletes any `Vec`s or `Mat`s allocated by the line search.

252:   Collective

254:   Input Parameter:
255: . linesearch - The `SNESLineSearch` instance.

257:   Level: developer

259:   Note:
260:   Usually only called by `SNESReset()`

262: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchSetUp()`
263: @*/
264: PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch)
265: {
266:   PetscFunctionBegin;
267:   if (linesearch->ops->reset) PetscUseTypeMethod(linesearch, reset);

269:   PetscCall(VecDestroy(&linesearch->vec_sol_new));
270:   PetscCall(VecDestroy(&linesearch->vec_func_new));

272:   PetscCall(VecDestroyVecs(linesearch->nwork, &linesearch->work));

274:   linesearch->nwork       = 0;
275:   linesearch->setupcalled = PETSC_FALSE;
276:   PetscFunctionReturn(PETSC_SUCCESS);
277: }

279: /*@C
280:   SNESLineSearchSetFunction - Sets the function evaluation used by the `SNES` line search
281:   `

283:   Input Parameters:
284: + linesearch - the `SNESLineSearch` context
285: - func       - function evaluation routine, this is usually the function provided with `SNESSetFunction()`

287:   Calling sequence of `func`:
288: + snes - the `SNES` with which the `SNESLineSearch` context is associated with
289: . x    - the input vector
290: - f    - the computed value of the function

292:   Level: developer

294:   Note:
295:   By default the `SNESLineSearch` uses the function provided by `SNESSetFunction()` so this is rarely needed

297: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESSetFunction()`
298: @*/
299: PetscErrorCode SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES snes, Vec x, Vec f))
300: {
301:   PetscFunctionBegin;
303:   linesearch->ops->snesfunc = func;
304:   PetscFunctionReturn(PETSC_SUCCESS);
305: }

307: /*@C
308:   SNESLineSearchSetPreCheck - Sets a function that is called after the initial search direction has been computed but
309:   before the line search routine has been applied. Allows adjusting the result of (usually a linear solve) that
310:   determined the search direction.

312:   Logically Collective

314:   Input Parameters:
315: + linesearch - the `SNESLineSearch` context
316: . func       - [optional] function evaluation routine
317: - ctx        - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)

319:   Calling sequence of `func`:
320: + ls        - the `SNESLineSearch` context
321: . x         - the current solution
322: . d         - the current search direction
323: . changed_d - indicates if the search direction has been changed
324: - ctx       - the context passed to `SNESLineSearchSetPreCheck()`

326:   Level: intermediate

328:   Note:
329:   Use `SNESLineSearchSetPostCheck()` to change the step after the line search is complete.

331:   Use `SNESVISetVariableBounds()` and `SNESVISetComputeVariableBounds()` to cause `SNES` to automatically control the ranges of variables allowed.

333: .seealso: [](ch_snes), `SNES`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()`,
334:           `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`

336: @*/
337: PetscErrorCode SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch ls, Vec x, Vec d, PetscBool *changed_d, void *ctx), void *ctx)
338: {
339:   PetscFunctionBegin;
341:   if (func) linesearch->ops->precheck = func;
342:   if (ctx) linesearch->precheckctx = ctx;
343:   PetscFunctionReturn(PETSC_SUCCESS);
344: }

346: /*@C
347:   SNESLineSearchGetPreCheck - Gets the pre-check function for the line search routine.

349:   Input Parameter:
350: . linesearch - the `SNESLineSearch` context

352:   Output Parameters:
353: + func - [optional] function evaluation routine,  for calling sequence see `SNESLineSearchSetPreCheck()`
354: - ctx  - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)

356:   Level: intermediate

358: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchSetPostCheck()`
359: @*/
360: PetscErrorCode SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch, Vec, Vec, PetscBool *, void *), void **ctx)
361: {
362:   PetscFunctionBegin;
364:   if (func) *func = linesearch->ops->precheck;
365:   if (ctx) *ctx = linesearch->precheckctx;
366:   PetscFunctionReturn(PETSC_SUCCESS);
367: }

369: /*@C
370:   SNESLineSearchSetPostCheck - Sets a user function that is called after the line search has been applied to determine the step
371:   direction and length. Allows the user a chance to change or override the decision of the line search routine

373:   Logically Collective

375:   Input Parameters:
376: + linesearch - the `SNESLineSearch` context
377: . func       - [optional] function evaluation routine
378: - ctx        - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)

380:   Calling sequence of `func`:
381: + ls        - the `SNESLineSearch` context
382: . x         - the current solution
383: . d         - the current search direction
384: . w         - $ w = x + lambda*d $ for some lambda
385: . changed_d - indicates if the search direction `d` has been changed
386: . changed_w - indicates `w` has been changed
387: - ctx       - the context passed to `SNESLineSearchSetPreCheck()`

389:   Level: intermediate

391:   Notes:
392:   Use `SNESLineSearchSetPreCheck()` to change the step before the line search is complete.

394:   Use `SNESVISetVariableBounds()` and `SNESVISetComputeVariableBounds()` to cause `SNES` to automatically control the ranges of variables allowed.

396: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchGetPostCheck()`,
397:           `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`
398: @*/
399: PetscErrorCode SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch ls, Vec x, Vec d, Vec w, PetscBool *changed_d, PetscBool *changed_w, void *ctx), void *ctx)
400: {
401:   PetscFunctionBegin;
403:   if (func) linesearch->ops->postcheck = func;
404:   if (ctx) linesearch->postcheckctx = ctx;
405:   PetscFunctionReturn(PETSC_SUCCESS);
406: }

408: /*@C
409:   SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine.

411:   Input Parameter:
412: . linesearch - the `SNESLineSearch` context

414:   Output Parameters:
415: + func - [optional] function evaluation routine, see for the calling sequence `SNESLineSearchSetPostCheck()`
416: - ctx  - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)

418:   Level: intermediate

420: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`
421: @*/
422: PetscErrorCode SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx)
423: {
424:   PetscFunctionBegin;
426:   if (func) *func = linesearch->ops->postcheck;
427:   if (ctx) *ctx = linesearch->postcheckctx;
428:   PetscFunctionReturn(PETSC_SUCCESS);
429: }

431: /*@
432:   SNESLineSearchPreCheck - Prepares the line search for being applied.

434:   Logically Collective

436:   Input Parameters:
437: + linesearch - The linesearch instance.
438: . X          - The current solution
439: - Y          - The step direction

441:   Output Parameter:
442: . changed - Indicator that the precheck routine has changed `Y`

444:   Level: advanced

446:   Note:
447:   This calls any function provided with `SNESLineSearchSetPreCheck()` and is called automatically inside the line search routines

449:   Developer Note:
450:   The use of `PetscObjectGetState()` would eliminate the need for the `changed` argument to be provided

452: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`,
453:           `SNESLineSearchGetPostCheck()`
454: @*/
455: PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed)
456: {
457:   PetscFunctionBegin;
458:   *changed = PETSC_FALSE;
459:   if (linesearch->ops->precheck) {
460:     PetscUseTypeMethod(linesearch, precheck, X, Y, changed, linesearch->precheckctx);
462:   }
463:   PetscFunctionReturn(PETSC_SUCCESS);
464: }

466: /*@
467:   SNESLineSearchPostCheck - Hook to modify step direction or updated solution after a successful linesearch

469:   Logically Collective

471:   Input Parameters:
472: + linesearch - The line search context
473: . X          - The last solution
474: . Y          - The step direction
475: - W          - The updated solution, $ W = X + lambda*Y $ for some lambda

477:   Output Parameters:
478: + changed_Y - Indicator if the direction `Y` has been changed.
479: - changed_W - Indicator if the new candidate solution `W` has been changed.

481:   Level: developer

483:   Note:
484:   This calls any function provided with `SNESLineSearchSetPostCheck()` and is called automatically inside the line search routines

486:   Developer Note:
487:   The use of `PetscObjectGetState()` would eliminate the need for the `changed_Y` and `changed_W` arguments to be provided

489: .seealso: [](ch_snes), `SNES`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPrecheck()`, `SNESLineSearchGetPrecheck()`
490: @*/
491: PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W)
492: {
493:   PetscFunctionBegin;
494:   *changed_Y = PETSC_FALSE;
495:   *changed_W = PETSC_FALSE;
496:   if (linesearch->ops->postcheck) {
497:     PetscUseTypeMethod(linesearch, postcheck, X, Y, W, changed_Y, changed_W, linesearch->postcheckctx);
500:   }
501:   PetscFunctionReturn(PETSC_SUCCESS);
502: }

504: /*@C
505:   SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration {cite}`hindmarsh1996time`

507:   Logically Collective

509:   Input Parameters:
510: + linesearch - the line search context
511: . X          - base state for this step
512: - ctx        - context for this function

514:   Input/Output Parameter:
515: . Y - correction, possibly modified

517:   Output Parameter:
518: . changed - flag indicating that `Y` was modified

520:   Options Database Keys:
521: + -snes_linesearch_precheck_picard       - activate this routine
522: - -snes_linesearch_precheck_picard_angle - angle

524:   Level: advanced

526:   Notes:
527:   This function should be passed to `SNESLineSearchSetPreCheck()`

529:   The justification for this method involves the linear convergence of a Picard iteration
530:   so the Picard linearization should be provided in place of the "Jacobian"  {cite}`hindmarsh1996time`. This correction
531:   is generally not useful when using a Newton linearization.

533:   Developer Note:
534:   The use of `PetscObjectGetState()` would eliminate the need for the `changed` argument to be provided

536: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESSetPicard()`, `SNESGetLineSearch()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchSetPostCheck()`
537: @*/
538: PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed, void *ctx)
539: {
540:   PetscReal   angle = *(PetscReal *)linesearch->precheckctx;
541:   Vec         Ylast;
542:   PetscScalar dot;
543:   PetscInt    iter;
544:   PetscReal   ynorm, ylastnorm, theta, angle_radians;
545:   SNES        snes;

547:   PetscFunctionBegin;
548:   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
549:   PetscCall(PetscObjectQuery((PetscObject)snes, "SNESLineSearchPreCheckPicard_Ylast", (PetscObject *)&Ylast));
550:   if (!Ylast) {
551:     PetscCall(VecDuplicate(Y, &Ylast));
552:     PetscCall(PetscObjectCompose((PetscObject)snes, "SNESLineSearchPreCheckPicard_Ylast", (PetscObject)Ylast));
553:     PetscCall(PetscObjectDereference((PetscObject)Ylast));
554:   }
555:   PetscCall(SNESGetIterationNumber(snes, &iter));
556:   if (iter < 2) {
557:     PetscCall(VecCopy(Y, Ylast));
558:     *changed = PETSC_FALSE;
559:     PetscFunctionReturn(PETSC_SUCCESS);
560:   }

562:   PetscCall(VecDot(Y, Ylast, &dot));
563:   PetscCall(VecNorm(Y, NORM_2, &ynorm));
564:   PetscCall(VecNorm(Ylast, NORM_2, &ylastnorm));
565:   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
566:   theta         = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm), -1.0, 1.0));
567:   angle_radians = angle * PETSC_PI / 180.;
568:   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
569:     /* Modify the step Y */
570:     PetscReal alpha, ydiffnorm;
571:     PetscCall(VecAXPY(Ylast, -1.0, Y));
572:     PetscCall(VecNorm(Ylast, NORM_2, &ydiffnorm));
573:     alpha = (ydiffnorm > .001 * ylastnorm) ? ylastnorm / ydiffnorm : 1000.0;
574:     PetscCall(VecCopy(Y, Ylast));
575:     PetscCall(VecScale(Y, alpha));
576:     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));
577:     *changed = PETSC_TRUE;
578:   } else {
579:     PetscCall(PetscInfo(snes, "Angle %14.12e degrees exceeds threshold %14.12e, no correction applied\n", (double)(theta * 180. / PETSC_PI), (double)angle));
580:     PetscCall(VecCopy(Y, Ylast));
581:     *changed = PETSC_FALSE;
582:   }
583:   PetscFunctionReturn(PETSC_SUCCESS);
584: }

586: /*@
587:   SNESLineSearchApply - Computes the line-search update.

589:   Collective

591:   Input Parameters:
592: + linesearch - The line search context
593: - Y          - The search direction

595:   Input/Output Parameters:
596: + X     - The current solution, on output the new solution
597: . F     - The current function value, on output the new function value at the solution value `X`
598: - fnorm - The current norm of `F`, on output the new norm of `F`

600:   Options Database Keys:
601: + -snes_linesearch_type                - basic (or equivalently none), bt, l2, cp, nleqerr, shell
602: . -snes_linesearch_monitor [:filename] - Print progress of line searches
603: . -snes_linesearch_damping             - The linesearch damping parameter, default is 1.0 (no damping)
604: . -snes_linesearch_norms               - Turn on/off the linesearch norms computation (SNESLineSearchSetComputeNorms())
605: . -snes_linesearch_keeplambda          - Keep the previous search length as the initial guess
606: - -snes_linesearch_max_it              - The number of iterations for iterative line searches

608:   Level: intermediate

610:   Notes:
611:   This is typically called from within a `SNESSolve()` implementation in order to
612:   help with convergence of the nonlinear method.  Various `SNES` types use line searches
613:   in different ways, but the overarching theme is that a line search is used to determine
614:   an optimal damping parameter of a step at each iteration of the method.  Each
615:   application of the line search may invoke `SNESComputeFunction()` several times, and
616:   therefore may be fairly expensive.

618: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchPreCheck()`, `SNESLineSearchPostCheck()`, `SNESSolve()`, `SNESComputeFunction()`, `SNESLineSearchSetComputeNorms()`,
619:           `SNESLineSearchType`, `SNESLineSearchSetType()`
620: @*/
621: PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal *fnorm, Vec Y)
622: {
623:   PetscFunctionBegin;

629:   linesearch->result = SNES_LINESEARCH_SUCCEEDED;

631:   linesearch->vec_sol    = X;
632:   linesearch->vec_update = Y;
633:   linesearch->vec_func   = F;

635:   PetscCall(SNESLineSearchSetUp(linesearch));

637:   if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */

639:   if (fnorm) linesearch->fnorm = *fnorm;
640:   else PetscCall(VecNorm(F, NORM_2, &linesearch->fnorm));

642:   PetscCall(PetscLogEventBegin(SNESLINESEARCH_Apply, linesearch, X, F, Y));

644:   PetscUseTypeMethod(linesearch, apply);

646:   PetscCall(PetscLogEventEnd(SNESLINESEARCH_Apply, linesearch, X, F, Y));

648:   if (fnorm) *fnorm = linesearch->fnorm;
649:   PetscFunctionReturn(PETSC_SUCCESS);
650: }

652: /*@
653:   SNESLineSearchDestroy - Destroys the line search instance.

655:   Collective

657:   Input Parameter:
658: . linesearch - The line search context

660:   Level: developer

662:   Note:
663:   The line search in `SNES` is automatically called on `SNESDestroy()` so this call is rarely needed

665: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchReset()`, `SNESDestroy()`
666: @*/
667: PetscErrorCode SNESLineSearchDestroy(SNESLineSearch *linesearch)
668: {
669:   PetscFunctionBegin;
670:   if (!*linesearch) PetscFunctionReturn(PETSC_SUCCESS);
672:   if (--((PetscObject)(*linesearch))->refct > 0) {
673:     *linesearch = NULL;
674:     PetscFunctionReturn(PETSC_SUCCESS);
675:   }
676:   PetscCall(PetscObjectSAWsViewOff((PetscObject)*linesearch));
677:   PetscCall(SNESLineSearchReset(*linesearch));
678:   PetscTryTypeMethod(*linesearch, destroy);
679:   PetscCall(PetscViewerDestroy(&(*linesearch)->monitor));
680:   PetscCall(SNESLineSearchMonitorCancel((*linesearch)));
681:   PetscCall(PetscHeaderDestroy(linesearch));
682:   PetscFunctionReturn(PETSC_SUCCESS);
683: }

685: /*@
686:   SNESLineSearchSetDefaultMonitor - Turns on/off printing useful information and debugging output about the line search.

688:   Logically Collective

690:   Input Parameters:
691: + linesearch - the linesearch object
692: - viewer     - an `PETSCVIEWERASCII` `PetscViewer` or `NULL` to turn off monitor

694:   Options Database Key:
695: . -snes_linesearch_monitor [:filename] - enables the monitor

697:   Level: intermediate

699:   Developer Notes:
700:   This monitor is implemented differently than the other line search monitors that are set with
701:   `SNESLineSearchMonitorSet()` since it is called in many locations of the line search routines to display aspects of the
702:   line search that are not visible to the other monitors.

704: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `PETSCVIEWERASCII`, `SNESGetLineSearch()`, `SNESLineSearchGetDefaultMonitor()`, `PetscViewer`, `SNESLineSearchSetMonitor()`,
705:           `SNESLineSearchMonitorSetFromOptions()`
706: @*/
707: PetscErrorCode SNESLineSearchSetDefaultMonitor(SNESLineSearch linesearch, PetscViewer viewer)
708: {
709:   PetscFunctionBegin;
710:   if (viewer) PetscCall(PetscObjectReference((PetscObject)viewer));
711:   PetscCall(PetscViewerDestroy(&linesearch->monitor));
712:   linesearch->monitor = viewer;
713:   PetscFunctionReturn(PETSC_SUCCESS);
714: }

716: /*@
717:   SNESLineSearchGetDefaultMonitor - Gets the `PetscViewer` instance for the default line search monitor that is turned on with `SNESLineSearchSetDefaultMonitor()`

719:   Logically Collective

721:   Input Parameter:
722: . linesearch - the line search context

724:   Output Parameter:
725: . monitor - monitor context

727:   Level: intermediate

729: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchSetDefaultMonitor()`, `PetscViewer`
730: @*/
731: PetscErrorCode SNESLineSearchGetDefaultMonitor(SNESLineSearch linesearch, PetscViewer *monitor)
732: {
733:   PetscFunctionBegin;
735:   *monitor = linesearch->monitor;
736:   PetscFunctionReturn(PETSC_SUCCESS);
737: }

739: /*@C
740:   SNESLineSearchMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated in the options database

742:   Collective

744:   Input Parameters:
745: + ls           - `SNESLineSearch` object to monitor
746: . name         - the monitor type
747: . help         - message indicating what monitoring is done
748: . manual       - manual page for the monitor
749: . monitor      - the monitor function
750: - monitorsetup - a function that is called once ONLY if the user selected this monitor that may set additional features of the `SNESLineSearch` or `PetscViewer`

752:   Calling sequence of `monitor`:
753: + ls - `SNESLineSearch` object being monitored
754: - vf - a `PetscViewerAndFormat` struct that provides the `PetscViewer` and `PetscViewerFormat` being used

756:   Calling sequence of `monitorsetup`:
757: + ls - `SNESLineSearch` object being monitored
758: - vf - a `PetscViewerAndFormat` struct that provides the `PetscViewer` and `PetscViewerFormat` being used

760:   Level: advanced

762: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetMonitor()`, `PetscOptionsGetViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
763:           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
764:           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
765:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
766:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
767:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
768:           `PetscOptionsFList()`, `PetscOptionsEList()`
769: @*/
770: PetscErrorCode SNESLineSearchMonitorSetFromOptions(SNESLineSearch ls, const char name[], const char help[], const char manual[], PetscErrorCode (*monitor)(SNESLineSearch ls, PetscViewerAndFormat *vf), PetscErrorCode (*monitorsetup)(SNESLineSearch ls, PetscViewerAndFormat *vf))
771: {
772:   PetscViewer       viewer;
773:   PetscViewerFormat format;
774:   PetscBool         flg;

776:   PetscFunctionBegin;
777:   PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)ls), ((PetscObject)ls)->options, ((PetscObject)ls)->prefix, name, &viewer, &format, &flg));
778:   if (flg) {
779:     PetscViewerAndFormat *vf;
780:     PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf));
781:     PetscCall(PetscObjectDereference((PetscObject)viewer));
782:     if (monitorsetup) PetscCall((*monitorsetup)(ls, vf));
783:     PetscCall(SNESLineSearchMonitorSet(ls, (PetscErrorCode(*)(SNESLineSearch, void *))monitor, vf, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy));
784:   }
785:   PetscFunctionReturn(PETSC_SUCCESS);
786: }

788: /*@
789:   SNESLineSearchSetFromOptions - Sets options for the line search

791:   Logically Collective

793:   Input Parameter:
794: . linesearch - a `SNESLineSearch` line search context

796:   Options Database Keys:
797: + -snes_linesearch_type <type>                                      - basic (or equivalently none), bt, l2, cp, nleqerr, shell
798: . -snes_linesearch_order <order>                                    - 1, 2, 3.  Most types only support certain orders (bt supports 2 or 3)
799: . -snes_linesearch_norms                                            - Turn on/off the linesearch norms for the basic linesearch typem (`SNESLineSearchSetComputeNorms()`)
800: . -snes_linesearch_minlambda                                        - The minimum step length
801: . -snes_linesearch_maxstep                                          - The maximum step size
802: . -snes_linesearch_rtol                                             - Relative tolerance for iterative line searches
803: . -snes_linesearch_atol                                             - Absolute tolerance for iterative line searches
804: . -snes_linesearch_ltol                                             - Change in lambda tolerance for iterative line searches
805: . -snes_linesearch_max_it                                           - The number of iterations for iterative line searches
806: . -snes_linesearch_monitor [:filename]                              - Print progress of line searches
807: . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine
808: . -snes_linesearch_damping                                          - The linesearch damping parameter
809: . -snes_linesearch_keeplambda                                       - Keep the previous search length as the initial guess.
810: . -snes_linesearch_precheck_picard                                  - Use precheck that speeds up convergence of picard method
811: - -snes_linesearch_precheck_picard_angle                            - Angle used in Picard precheck method

813:   Level: intermediate

815: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchCreate()`, `SNESLineSearchSetOrder()`, `SNESLineSearchSetType()`, `SNESLineSearchSetTolerances()`, `SNESLineSearchSetDamping()`, `SNESLineSearchPreCheckPicard()`,
816:           `SNESLineSearchType`, `SNESLineSearchSetComputeNorms()`
817: @*/
818: PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch)
819: {
820:   const char *deft = SNESLINESEARCHBASIC;
821:   char        type[256];
822:   PetscBool   flg, set;
823:   PetscViewer viewer;

825:   PetscFunctionBegin;
826:   PetscCall(SNESLineSearchRegisterAll());

828:   PetscObjectOptionsBegin((PetscObject)linesearch);
829:   if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name;
830:   PetscCall(PetscOptionsFList("-snes_linesearch_type", "Linesearch type", "SNESLineSearchSetType", SNESLineSearchList, deft, type, 256, &flg));
831:   if (flg) {
832:     PetscCall(SNESLineSearchSetType(linesearch, type));
833:   } else if (!((PetscObject)linesearch)->type_name) {
834:     PetscCall(SNESLineSearchSetType(linesearch, deft));
835:   }

837:   PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)linesearch), ((PetscObject)linesearch)->options, ((PetscObject)linesearch)->prefix, "-snes_linesearch_monitor", &viewer, NULL, &set));
838:   if (set) {
839:     PetscCall(SNESLineSearchSetDefaultMonitor(linesearch, viewer));
840:     PetscCall(PetscViewerDestroy(&viewer));
841:   }
842:   PetscCall(SNESLineSearchMonitorSetFromOptions(linesearch, "-snes_linesearch_monitor_solution_update", "View correction at each iteration", "SNESLineSearchMonitorSolutionUpdate", SNESLineSearchMonitorSolutionUpdate, NULL));

844:   /* tolerances */
845:   PetscCall(PetscOptionsReal("-snes_linesearch_minlambda", "Minimum step length", "SNESLineSearchSetTolerances", linesearch->steptol, &linesearch->steptol, NULL));
846:   PetscCall(PetscOptionsReal("-snes_linesearch_maxstep", "Maximum step size", "SNESLineSearchSetTolerances", linesearch->maxstep, &linesearch->maxstep, NULL));
847:   PetscCall(PetscOptionsReal("-snes_linesearch_rtol", "Relative tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->rtol, &linesearch->rtol, NULL));
848:   PetscCall(PetscOptionsReal("-snes_linesearch_atol", "Absolute tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->atol, &linesearch->atol, NULL));
849:   PetscCall(PetscOptionsReal("-snes_linesearch_ltol", "Change in lambda tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->ltol, &linesearch->ltol, NULL));
850:   PetscCall(PetscOptionsInt("-snes_linesearch_max_it", "Maximum iterations for iterative line searches", "SNESLineSearchSetTolerances", linesearch->max_its, &linesearch->max_its, NULL));

852:   /* damping parameters */
853:   PetscCall(PetscOptionsReal("-snes_linesearch_damping", "Line search damping and initial step guess", "SNESLineSearchSetDamping", linesearch->damping, &linesearch->damping, NULL));

855:   PetscCall(PetscOptionsBool("-snes_linesearch_keeplambda", "Use previous lambda as damping", "SNESLineSearchSetKeepLambda", linesearch->keeplambda, &linesearch->keeplambda, NULL));

857:   /* precheck */
858:   PetscCall(PetscOptionsBool("-snes_linesearch_precheck_picard", "Use a correction that sometimes improves convergence of Picard iteration", "SNESLineSearchPreCheckPicard", flg, &flg, &set));
859:   if (set) {
860:     if (flg) {
861:       linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */

863:       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));
864:       PetscCall(SNESLineSearchSetPreCheck(linesearch, SNESLineSearchPreCheckPicard, &linesearch->precheck_picard_angle));
865:     } else {
866:       PetscCall(SNESLineSearchSetPreCheck(linesearch, NULL, NULL));
867:     }
868:   }
869:   PetscCall(PetscOptionsInt("-snes_linesearch_order", "Order of approximation used in the line search", "SNESLineSearchSetOrder", linesearch->order, &linesearch->order, NULL));
870:   PetscCall(PetscOptionsBool("-snes_linesearch_norms", "Compute final norms in line search", "SNESLineSearchSetComputeNorms", linesearch->norms, &linesearch->norms, NULL));

872:   PetscTryTypeMethod(linesearch, setfromoptions, PetscOptionsObject);

874:   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)linesearch, PetscOptionsObject));
875:   PetscOptionsEnd();
876:   PetscFunctionReturn(PETSC_SUCCESS);
877: }

879: /*@
880:   SNESLineSearchView - Prints useful information about the line search

882:   Logically Collective

884:   Input Parameters:
885: + linesearch - line search context
886: - viewer     - the `PetscViewer` to display the line search information to

888:   Level: intermediate

890: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `PetscViewer`, `SNESLineSearchCreate()`
891: @*/
892: PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer)
893: {
894:   PetscBool iascii;

896:   PetscFunctionBegin;
898:   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch), &viewer));
900:   PetscCheckSameComm(linesearch, 1, viewer, 2);

902:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
903:   if (iascii) {
904:     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)linesearch, viewer));
905:     PetscCall(PetscViewerASCIIPushTab(viewer));
906:     PetscTryTypeMethod(linesearch, view, viewer);
907:     PetscCall(PetscViewerASCIIPopTab(viewer));
908:     PetscCall(PetscViewerASCIIPrintf(viewer, "  maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep, (double)linesearch->steptol));
909:     PetscCall(PetscViewerASCIIPrintf(viewer, "  tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol, (double)linesearch->atol, (double)linesearch->ltol));
910:     PetscCall(PetscViewerASCIIPrintf(viewer, "  maximum iterations=%" PetscInt_FMT "\n", linesearch->max_its));
911:     if (linesearch->ops->precheck) {
912:       if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) {
913:         PetscCall(PetscViewerASCIIPrintf(viewer, "  using precheck step to speed up Picard convergence\n"));
914:       } else {
915:         PetscCall(PetscViewerASCIIPrintf(viewer, "  using user-defined precheck step\n"));
916:       }
917:     }
918:     if (linesearch->ops->postcheck) PetscCall(PetscViewerASCIIPrintf(viewer, "  using user-defined postcheck step\n"));
919:   }
920:   PetscFunctionReturn(PETSC_SUCCESS);
921: }

923: /*@C
924:   SNESLineSearchGetType - Gets the `SNESLinesearchType` of a `SNESLineSearch`

926:   Logically Collective

928:   Input Parameter:
929: . linesearch - the line search context

931:   Output Parameter:
932: . type - The type of line search, or `NULL` if not set

934:   Level: intermediate

936: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchSetType()`
937: @*/
938: PetscErrorCode SNESLineSearchGetType(SNESLineSearch linesearch, SNESLineSearchType *type)
939: {
940:   PetscFunctionBegin;
942:   PetscAssertPointer(type, 2);
943:   *type = ((PetscObject)linesearch)->type_name;
944:   PetscFunctionReturn(PETSC_SUCCESS);
945: }

947: /*@C
948:   SNESLineSearchSetType - Sets the `SNESLinesearchType` of a `SNESLineSearch`

950:   Logically Collective

952:   Input Parameters:
953: + linesearch - the line search context
954: - type       - The type of line search to be used, see `SNESLineSearchType`

956:   Options Database Key:
957: . -snes_linesearch_type <type> - basic (or equivalently none), bt, l2, cp, nleqerr, shell

959:   Level: intermediate

961: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchGetType()`
962: @*/
963: PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type)
964: {
965:   PetscBool match;
966:   PetscErrorCode (*r)(SNESLineSearch);

968:   PetscFunctionBegin;
970:   PetscAssertPointer(type, 2);

972:   PetscCall(PetscObjectTypeCompare((PetscObject)linesearch, type, &match));
973:   if (match) PetscFunctionReturn(PETSC_SUCCESS);

975:   PetscCall(PetscFunctionListFind(SNESLineSearchList, type, &r));
976:   PetscCheck(r, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested Line Search type %s", type);
977:   /* Destroy the previous private line search context */
978:   if (linesearch->ops->destroy) {
979:     PetscCall((*(linesearch)->ops->destroy)(linesearch));
980:     linesearch->ops->destroy = NULL;
981:   }
982:   /* Reinitialize function pointers in SNESLineSearchOps structure */
983:   linesearch->ops->apply          = NULL;
984:   linesearch->ops->view           = NULL;
985:   linesearch->ops->setfromoptions = NULL;
986:   linesearch->ops->destroy        = NULL;

988:   PetscCall(PetscObjectChangeTypeName((PetscObject)linesearch, type));
989:   PetscCall((*r)(linesearch));
990:   PetscFunctionReturn(PETSC_SUCCESS);
991: }

993: /*@
994:   SNESLineSearchSetSNES - Sets the `SNES` for the linesearch for function evaluation.

996:   Input Parameters:
997: + linesearch - the line search context
998: - snes       - The `SNES` instance

1000:   Level: developer

1002:   Note:
1003:   This happens automatically when the line search is obtained/created with
1004:   `SNESGetLineSearch()`.  This routine is therefore mainly called within `SNES`
1005:   implementations.

1007: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetSNES()`, `SNESLineSearchSetVecs()`
1008: @*/
1009: PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes)
1010: {
1011:   PetscFunctionBegin;
1014:   linesearch->snes = snes;
1015:   PetscFunctionReturn(PETSC_SUCCESS);
1016: }

1018: /*@
1019:   SNESLineSearchGetSNES - Gets the `SNES` instance associated with the line search.

1021:   Not Collective

1023:   Input Parameter:
1024: . linesearch - the line search context

1026:   Output Parameter:
1027: . snes - The `SNES` instance

1029:   Level: developer

1031: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESType`, `SNESLineSearchSetVecs()`
1032: @*/
1033: PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes)
1034: {
1035:   PetscFunctionBegin;
1037:   PetscAssertPointer(snes, 2);
1038:   *snes = linesearch->snes;
1039:   PetscFunctionReturn(PETSC_SUCCESS);
1040: }

1042: /*@
1043:   SNESLineSearchGetLambda - Gets the last line search steplength used

1045:   Not Collective

1047:   Input Parameter:
1048: . linesearch - the line search context

1050:   Output Parameter:
1051: . lambda - The last steplength computed during `SNESLineSearchApply()`

1053:   Level: advanced

1055:   Note:
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: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetLambda()`, `SNESLineSearchGetDamping()`, `SNESLineSearchApply()`
1062: @*/
1063: PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch linesearch, PetscReal *lambda)
1064: {
1065:   PetscFunctionBegin;
1067:   PetscAssertPointer(lambda, 2);
1068:   *lambda = linesearch->lambda;
1069:   PetscFunctionReturn(PETSC_SUCCESS);
1070: }

1072: /*@
1073:   SNESLineSearchSetLambda - Sets the line search steplength

1075:   Input Parameters:
1076: + linesearch - line search context
1077: - lambda     - The steplength to use

1079:   Level: advanced

1081:   Note:
1082:   This routine is typically used within implementations of `SNESLineSearchApply()`
1083:   to set the final steplength.  This routine (and `SNESLineSearchGetLambda()`) were
1084:   added in order to facilitate Quasi-Newton methods that use the previous steplength
1085:   as an inner scaling parameter.

1087: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetLambda()`
1088: @*/
1089: PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
1090: {
1091:   PetscFunctionBegin;
1093:   linesearch->lambda = lambda;
1094:   PetscFunctionReturn(PETSC_SUCCESS);
1095: }

1097: /*@
1098:   SNESLineSearchGetTolerances - Gets the tolerances for the line search.

1100:   Not Collective

1102:   Input Parameter:
1103: . linesearch - the line search context

1105:   Output Parameters:
1106: + steptol - The minimum steplength
1107: . maxstep - The maximum steplength
1108: . rtol    - The relative tolerance for iterative line searches
1109: . atol    - The absolute tolerance for iterative line searches
1110: . ltol    - The change in lambda tolerance for iterative line searches
1111: - max_its - The maximum number of iterations of the line search

1113:   Level: intermediate

1115:   Note:
1116:   Different line searches may implement these parameters slightly differently as
1117:   the type requires.

1119: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetTolerances()`
1120: @*/
1121: PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch linesearch, PetscReal *steptol, PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
1122: {
1123:   PetscFunctionBegin;
1125:   if (steptol) {
1126:     PetscAssertPointer(steptol, 2);
1127:     *steptol = linesearch->steptol;
1128:   }
1129:   if (maxstep) {
1130:     PetscAssertPointer(maxstep, 3);
1131:     *maxstep = linesearch->maxstep;
1132:   }
1133:   if (rtol) {
1134:     PetscAssertPointer(rtol, 4);
1135:     *rtol = linesearch->rtol;
1136:   }
1137:   if (atol) {
1138:     PetscAssertPointer(atol, 5);
1139:     *atol = linesearch->atol;
1140:   }
1141:   if (ltol) {
1142:     PetscAssertPointer(ltol, 6);
1143:     *ltol = linesearch->ltol;
1144:   }
1145:   if (max_its) {
1146:     PetscAssertPointer(max_its, 7);
1147:     *max_its = linesearch->max_its;
1148:   }
1149:   PetscFunctionReturn(PETSC_SUCCESS);
1150: }

1152: /*@
1153:   SNESLineSearchSetTolerances -  Sets the tolerances for the linesearch.

1155:   Collective

1157:   Input Parameters:
1158: + linesearch - the line search context
1159: . steptol    - The minimum steplength
1160: . maxstep    - The maximum steplength
1161: . rtol       - The relative tolerance for iterative line searches
1162: . atol       - The absolute tolerance for iterative line searches
1163: . ltol       - The change in lambda tolerance for iterative line searches
1164: - max_it     - The maximum number of iterations of the line search

1166:   Options Database Keys:
1167: + -snes_linesearch_minlambda - The minimum step length
1168: . -snes_linesearch_maxstep   - The maximum step size
1169: . -snes_linesearch_rtol      - Relative tolerance for iterative line searches
1170: . -snes_linesearch_atol      - Absolute tolerance for iterative line searches
1171: . -snes_linesearch_ltol      - Change in lambda tolerance for iterative line searches
1172: - -snes_linesearch_max_it    - The number of iterations for iterative line searches

1174:   Level: intermediate

1176:   Note:
1177:   The user may choose to not set any of the tolerances using `PETSC_DEFAULT` in place of an argument.

1179: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetTolerances()`
1180: @*/
1181: PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch linesearch, PetscReal steptol, PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_it)
1182: {
1183:   PetscFunctionBegin;

1192:   if (steptol != (PetscReal)PETSC_DEFAULT) {
1193:     PetscCheck(steptol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Minimum step length %14.12e must be non-negative", (double)steptol);
1194:     linesearch->steptol = steptol;
1195:   }

1197:   if (maxstep != (PetscReal)PETSC_DEFAULT) {
1198:     PetscCheck(maxstep >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum step length %14.12e must be non-negative", (double)maxstep);
1199:     linesearch->maxstep = maxstep;
1200:   }

1202:   if (rtol != (PetscReal)PETSC_DEFAULT) {
1203:     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);
1204:     linesearch->rtol = rtol;
1205:   }

1207:   if (atol != (PetscReal)PETSC_DEFAULT) {
1208:     PetscCheck(atol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Absolute tolerance %14.12e must be non-negative", (double)atol);
1209:     linesearch->atol = atol;
1210:   }

1212:   if (ltol != (PetscReal)PETSC_DEFAULT) {
1213:     PetscCheck(ltol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Lambda tolerance %14.12e must be non-negative", (double)ltol);
1214:     linesearch->ltol = ltol;
1215:   }

1217:   if (max_it != PETSC_DEFAULT) {
1218:     PetscCheck(max_it >= 0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of iterations %" PetscInt_FMT " must be non-negative", max_it);
1219:     linesearch->max_its = max_it;
1220:   }
1221:   PetscFunctionReturn(PETSC_SUCCESS);
1222: }

1224: /*@
1225:   SNESLineSearchGetDamping - Gets the line search damping parameter.

1227:   Input Parameter:
1228: . linesearch - the line search context

1230:   Output Parameter:
1231: . damping - The damping parameter

1233:   Level: advanced

1235: .seealso: [](ch_snes), `SNES`, `SNESLineSearchGetStepTolerance()`, `SNESQN`
1236: @*/
1237: PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch linesearch, PetscReal *damping)
1238: {
1239:   PetscFunctionBegin;
1241:   PetscAssertPointer(damping, 2);
1242:   *damping = linesearch->damping;
1243:   PetscFunctionReturn(PETSC_SUCCESS);
1244: }

1246: /*@
1247:   SNESLineSearchSetDamping - Sets the line search damping parameter.

1249:   Input Parameters:
1250: + linesearch - the line search context
1251: - damping    - The damping parameter

1253:   Options Database Key:
1254: . -snes_linesearch_damping <damping> - the damping value

1256:   Level: intermediate

1258:   Note:
1259:   The `SNESLINESEARCHNONE` line search merely takes the update step scaled by the damping parameter.
1260:   The use of the damping parameter in the `SNESLINESEARCHL2` and `SNESLINESEARCHCP` line searches is much more subtle;
1261:   it is used as a starting point in calculating the secant step. However, the eventual
1262:   step may be of greater length than the damping parameter.  In the `SNESLINESEARCHBT` line search it is
1263:   used as the maximum possible step length, as the `SNESLINESEARCHBT` line search only backtracks.

1265: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetDamping()`
1266: @*/
1267: PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch linesearch, PetscReal damping)
1268: {
1269:   PetscFunctionBegin;
1271:   linesearch->damping = damping;
1272:   PetscFunctionReturn(PETSC_SUCCESS);
1273: }

1275: /*@
1276:   SNESLineSearchGetOrder - Gets the line search approximation order.

1278:   Input Parameter:
1279: . linesearch - the line search context

1281:   Output Parameter:
1282: . order - The order

1284:   Level: intermediate

1286: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetOrder()`
1287: @*/
1288: PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch linesearch, PetscInt *order)
1289: {
1290:   PetscFunctionBegin;
1292:   PetscAssertPointer(order, 2);
1293:   *order = linesearch->order;
1294:   PetscFunctionReturn(PETSC_SUCCESS);
1295: }

1297: /*@
1298:   SNESLineSearchSetOrder - Sets the maximum order of the polynomial fit used in the line search

1300:   Input Parameters:
1301: + linesearch - the line search context
1302: - order      - The order

1304:   Level: intermediate

1306:   Values for `order`\:
1307: +  1 or `SNES_LINESEARCH_ORDER_LINEAR` - linear order
1308: .  2 or `SNES_LINESEARCH_ORDER_QUADRATIC` - quadratic order
1309: -  3 or `SNES_LINESEARCH_ORDER_CUBIC` - cubic order

1311:   Options Database Key:
1312: . -snes_linesearch_order <order> - 1, 2, 3.  Most types only support certain orders (`SNESLINESEARCHBT` supports 2 or 3)

1314:   Note:
1315:   These orders are supported by `SNESLINESEARCHBT` and `SNESLINESEARCHCP`

1317: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetOrder()`, `SNESLineSearchSetDamping()`
1318: @*/
1319: PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch linesearch, PetscInt order)
1320: {
1321:   PetscFunctionBegin;
1323:   linesearch->order = order;
1324:   PetscFunctionReturn(PETSC_SUCCESS);
1325: }

1327: /*@
1328:   SNESLineSearchGetNorms - Gets the norms for the current solution `X`, the current update `Y`, and the current function value `F`.

1330:   Not Collective

1332:   Input Parameter:
1333: . linesearch - the line search context

1335:   Output Parameters:
1336: + xnorm - The norm of the current solution
1337: . fnorm - The norm of the current function, this is the `norm(function(X))` where `X` is the current solution.
1338: - ynorm - The norm of the current update

1340:   Level: developer

1342:   Notes:
1343:   Some values may not be up-to-date at particular points in the code.

1345:   This, in combination with `SNESLineSearchSetNorms()`, allow the line search and the `SNESSolve_XXX()` to share
1346:   computed values.

1348: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetNorms()` `SNESLineSearchGetVecs()`
1349: @*/
1350: PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal *xnorm, PetscReal *fnorm, PetscReal *ynorm)
1351: {
1352:   PetscFunctionBegin;
1354:   if (xnorm) *xnorm = linesearch->xnorm;
1355:   if (fnorm) *fnorm = linesearch->fnorm;
1356:   if (ynorm) *ynorm = linesearch->ynorm;
1357:   PetscFunctionReturn(PETSC_SUCCESS);
1358: }

1360: /*@
1361:   SNESLineSearchSetNorms - Gets the computed norms for the current solution `X`, the current update `Y`, and the current function value `F`.

1363:   Collective

1365:   Input Parameters:
1366: + linesearch - the line search context
1367: . xnorm      - The norm of the current solution
1368: . fnorm      - The norm of the current function, this is the `norm(function(X))` where `X` is the current solution
1369: - ynorm      - The norm of the current update

1371:   Level: developer

1373:   Note:
1374:   This is called by the line search routines to store the values they have just computed

1376: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()`
1377: @*/
1378: PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
1379: {
1380:   PetscFunctionBegin;
1382:   linesearch->xnorm = xnorm;
1383:   linesearch->fnorm = fnorm;
1384:   linesearch->ynorm = ynorm;
1385:   PetscFunctionReturn(PETSC_SUCCESS);
1386: }

1388: /*@
1389:   SNESLineSearchComputeNorms - Explicitly computes the norms of the current solution `X`, the current update `Y`, and the current function value `F`.

1391:   Input Parameter:
1392: . linesearch - the line search context

1394:   Options Database Key:
1395: . -snes_linesearch_norms - turn norm computation on or off

1397:   Level: intermediate

1399:   Developer Note:
1400:   The options database key is misnamed. It should be -snes_linesearch_compute_norms

1402: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms`, `SNESLineSearchSetNorms()`, `SNESLineSearchSetComputeNorms()`
1403: @*/
1404: PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
1405: {
1406:   SNES snes;

1408:   PetscFunctionBegin;
1409:   if (linesearch->norms) {
1410:     if (linesearch->ops->vinorm) {
1411:       PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
1412:       PetscCall(VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
1413:       PetscCall(VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm));
1414:       PetscCall((*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm));
1415:     } else {
1416:       PetscCall(VecNormBegin(linesearch->vec_func, NORM_2, &linesearch->fnorm));
1417:       PetscCall(VecNormBegin(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
1418:       PetscCall(VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm));
1419:       PetscCall(VecNormEnd(linesearch->vec_func, NORM_2, &linesearch->fnorm));
1420:       PetscCall(VecNormEnd(linesearch->vec_sol, NORM_2, &linesearch->xnorm));
1421:       PetscCall(VecNormEnd(linesearch->vec_update, NORM_2, &linesearch->ynorm));
1422:     }
1423:   }
1424:   PetscFunctionReturn(PETSC_SUCCESS);
1425: }

1427: /*@
1428:   SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search.

1430:   Input Parameters:
1431: + linesearch - the line search context
1432: - flg        - indicates whether or not to compute norms

1434:   Options Database Key:
1435: . -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic (none) `SNESLINESEARCHBASIC` line search

1437:   Level: intermediate

1439:   Note:
1440:   This is most relevant to the `SNESLINESEARCHBASIC` (or equivalently `SNESLINESEARCHNONE`) line search type since most line searches have a stopping criteria involving the norm.

1442:   Developer Note:
1443:   The options database key is misnamed. It should be -snes_linesearch_compute_norms

1445: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetNorms()`, `SNESLineSearchComputeNorms()`, `SNESLINESEARCHBASIC`
1446: @*/
1447: PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg)
1448: {
1449:   PetscFunctionBegin;
1450:   linesearch->norms = flg;
1451:   PetscFunctionReturn(PETSC_SUCCESS);
1452: }

1454: /*@
1455:   SNESLineSearchGetVecs - Gets the vectors from the `SNESLineSearch` context

1457:   Not Collective but the vectors are parallel

1459:   Input Parameter:
1460: . linesearch - the line search context

1462:   Output Parameters:
1463: + X - Solution vector
1464: . F - Function vector
1465: . Y - Search direction vector
1466: . W - Solution work vector
1467: - G - Function work vector

1469:   Level: advanced

1471:   Notes:
1472:   At the beginning of a line search application, `X` should contain a
1473:   solution and the vector `F` the function computed at `X`.  At the end of the
1474:   line search application, `X` should contain the new solution, and `F` the
1475:   function evaluated at the new solution.

1477:   These vectors are owned by the `SNESLineSearch` and should not be destroyed by the caller

1479: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()`
1480: @*/
1481: PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch, Vec *X, Vec *F, Vec *Y, Vec *W, Vec *G)
1482: {
1483:   PetscFunctionBegin;
1485:   if (X) {
1486:     PetscAssertPointer(X, 2);
1487:     *X = linesearch->vec_sol;
1488:   }
1489:   if (F) {
1490:     PetscAssertPointer(F, 3);
1491:     *F = linesearch->vec_func;
1492:   }
1493:   if (Y) {
1494:     PetscAssertPointer(Y, 4);
1495:     *Y = linesearch->vec_update;
1496:   }
1497:   if (W) {
1498:     PetscAssertPointer(W, 5);
1499:     *W = linesearch->vec_sol_new;
1500:   }
1501:   if (G) {
1502:     PetscAssertPointer(G, 6);
1503:     *G = linesearch->vec_func_new;
1504:   }
1505:   PetscFunctionReturn(PETSC_SUCCESS);
1506: }

1508: /*@
1509:   SNESLineSearchSetVecs - Sets the vectors on the `SNESLineSearch` context

1511:   Logically Collective

1513:   Input Parameters:
1514: + linesearch - the line search context
1515: . X          - Solution vector
1516: . F          - Function vector
1517: . Y          - Search direction vector
1518: . W          - Solution work vector
1519: - G          - Function work vector

1521:   Level: developer

1523: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetNorms()`, `SNESLineSearchGetVecs()`
1524: @*/
1525: PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch, Vec X, Vec F, Vec Y, Vec W, Vec G)
1526: {
1527:   PetscFunctionBegin;
1529:   if (X) {
1531:     linesearch->vec_sol = X;
1532:   }
1533:   if (F) {
1535:     linesearch->vec_func = F;
1536:   }
1537:   if (Y) {
1539:     linesearch->vec_update = Y;
1540:   }
1541:   if (W) {
1543:     linesearch->vec_sol_new = W;
1544:   }
1545:   if (G) {
1547:     linesearch->vec_func_new = G;
1548:   }
1549:   PetscFunctionReturn(PETSC_SUCCESS);
1550: }

1552: /*@C
1553:   SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1554:   `SNESLineSearch` options in the database.

1556:   Logically Collective

1558:   Input Parameters:
1559: + linesearch - the `SNESLineSearch` context
1560: - prefix     - the prefix to prepend to all option names

1562:   Level: advanced

1564:   Note:
1565:   A hyphen (-) must NOT be given at the beginning of the prefix name.
1566:   The first character of all runtime options is AUTOMATICALLY the hyphen.

1568: .seealso: [](ch_snes), `SNES`, `SNESLineSearch()`, `SNESLineSearchSetFromOptions()`, `SNESGetOptionsPrefix()`
1569: @*/
1570: PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch, const char prefix[])
1571: {
1572:   PetscFunctionBegin;
1574:   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)linesearch, prefix));
1575:   PetscFunctionReturn(PETSC_SUCCESS);
1576: }

1578: /*@C
1579:   SNESLineSearchGetOptionsPrefix - Gets the prefix used for searching for all
1580:   SNESLineSearch options in the database.

1582:   Not Collective

1584:   Input Parameter:
1585: . linesearch - the `SNESLineSearch` context

1587:   Output Parameter:
1588: . prefix - pointer to the prefix string used

1590:   Level: advanced

1592:   Fortran Notes:
1593:   The user should pass in a string 'prefix' of
1594:   sufficient length to hold the prefix.

1596: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESAppendOptionsPrefix()`
1597: @*/
1598: PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch, const char *prefix[])
1599: {
1600:   PetscFunctionBegin;
1602:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)linesearch, prefix));
1603:   PetscFunctionReturn(PETSC_SUCCESS);
1604: }

1606: /*@C
1607:   SNESLineSearchSetWorkVecs - Sets work vectors for the line search.

1609:   Input Parameters:
1610: + linesearch - the `SNESLineSearch` context
1611: - nwork      - the number of work vectors

1613:   Level: developer

1615:   Developer Note:
1616:   This is called from within the set up routines for each of the line search types `SNESLineSearchType`

1618: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESSetWorkVecs()`
1619: @*/
1620: PetscErrorCode SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork)
1621: {
1622:   PetscFunctionBegin;
1623:   PetscCheck(linesearch->vec_sol, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
1624:   PetscCall(VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work));
1625:   PetscFunctionReturn(PETSC_SUCCESS);
1626: }

1628: /*@
1629:   SNESLineSearchGetReason - Gets the success/failure status of the last line search application

1631:   Input Parameter:
1632: . linesearch - the line search context

1634:   Output Parameter:
1635: . result - The success or failure status

1637:   Level: developer

1639:   Note:
1640:   This is typically called after `SNESLineSearchApply()` in order to determine if the line search failed
1641:   (and set into the `SNES` convergence accordingly).

1643: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetReason()`, `SNESLineSearchReason`
1644: @*/
1645: PetscErrorCode SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result)
1646: {
1647:   PetscFunctionBegin;
1649:   PetscAssertPointer(result, 2);
1650:   *result = linesearch->result;
1651:   PetscFunctionReturn(PETSC_SUCCESS);
1652: }

1654: /*@
1655:   SNESLineSearchSetReason - Sets the success/failure status of the line search application

1657:   Input Parameters:
1658: + linesearch - the line search context
1659: - result     - The success or failure status

1661:   Level: developer

1663:   Note:
1664:   This is typically called in a `SNESLineSearchType` implementation of `SNESLineSearchApply()` or a `SNESLINESEARCHSHELL` implementation to set
1665:   the success or failure of the line search method.

1667: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchReason`, `SNESLineSearchGetSResult()`
1668: @*/
1669: PetscErrorCode SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result)
1670: {
1671:   PetscFunctionBegin;
1673:   linesearch->result = result;
1674:   PetscFunctionReturn(PETSC_SUCCESS);
1675: }

1677: // PetscClangLinter pragma disable: -fdoc-param-list-func-parameter-documentation
1678: /*@C
1679:   SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.

1681:   Logically Collective

1683:   Input Parameters:
1684: + linesearch  - the linesearch object
1685: . projectfunc - function for projecting the function to the bounds, see `SNESLineSearchSetVIFunctions` for calling sequence
1686: - normfunc    - function for computing the norm of an active set, see `SNESLineSearchSetVIFunctions ` for calling sequence

1688:   Level: advanced

1690:   Notes:
1691:   The VI solvers require projection of the solution to the feasible set.  `projectfunc` should implement this.

1693:   The VI solvers require special evaluation of the function norm such that the norm is only calculated
1694:   on the inactive set.  This should be implemented by `normfunc`.

1696: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetVIFunctions()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()`
1697: @*/
1698: PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc)
1699: {
1700:   PetscFunctionBegin;
1702:   if (projectfunc) linesearch->ops->viproject = projectfunc;
1703:   if (normfunc) linesearch->ops->vinorm = normfunc;
1704:   PetscFunctionReturn(PETSC_SUCCESS);
1705: }

1707: /*@C
1708:   SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.

1710:   Not Collective

1712:   Input Parameter:
1713: . linesearch - the line search context, obtain with `SNESGetLineSearch()`

1715:   Output Parameters:
1716: + projectfunc - function for projecting the function to the bounds
1717: - normfunc    - function for computing the norm of an active set

1719:   Level: advanced

1721: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetVIFunctions()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()`
1722: @*/
1723: PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc)
1724: {
1725:   PetscFunctionBegin;
1726:   if (projectfunc) *projectfunc = linesearch->ops->viproject;
1727:   if (normfunc) *normfunc = linesearch->ops->vinorm;
1728:   PetscFunctionReturn(PETSC_SUCCESS);
1729: }

1731: /*@C
1732:   SNESLineSearchRegister - register a line search type `SNESLineSearchType`

1734:   Input Parameters:
1735: + sname    - name of the `SNESLineSearchType()`
1736: - function - the creation function for that type

1738:   Calling sequence of `function`:
1739: . ls - the line search context

1741:   Level: advanced

1743: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchSetType()`
1744: @*/
1745: PetscErrorCode SNESLineSearchRegister(const char sname[], PetscErrorCode (*function)(SNESLineSearch ls))
1746: {
1747:   PetscFunctionBegin;
1748:   PetscCall(SNESInitializePackage());
1749:   PetscCall(PetscFunctionListAdd(&SNESLineSearchList, sname, function));
1750:   PetscFunctionReturn(PETSC_SUCCESS);
1751: }