Actual source code: linesearch.c

petsc-3.9.4 2018-09-11
Report Typos and Errors
  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: Only a single monitor function can be set for each SNESLineSearch object

100:    Level: intermediate

102: .keywords: SNESLineSearch, nonlinear, set, monitor

104: .seealso: SNESGetLineSearch(), SNESLineSearchMonitorDefault(), SNESLineSearchMonitorCancel()
105: @*/
106: PetscErrorCode  SNESLineSearchMonitorSet(SNESLineSearch ls,PetscErrorCode (*f)(SNESLineSearch,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
107: {
109:   PetscInt       i;
110:   PetscBool      identical;

114:   for (i=0; i<ls->numbermonitors;i++) {
115:     PetscMonitorCompare((PetscErrorCode (*)(void))f,mctx,monitordestroy,(PetscErrorCode (*)(void))ls->monitorftns[i],ls->monitorcontext[i],ls->monitordestroy[i],&identical);
116:     if (identical) return(0);
117:   }
118:   if (ls->numbermonitors >= MAXSNESLSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
119:   ls->monitorftns[ls->numbermonitors]          = f;
120:   ls->monitordestroy[ls->numbermonitors]   = monitordestroy;
121:   ls->monitorcontext[ls->numbermonitors++] = (void*)mctx;
122:   return(0);
123: }

125: /*@C
126:    SNESLineSearchMonitorSolutionUpdate - Monitors each update a new function value the linesearch tries

128:    Collective on SNESLineSearch

130:    Input Parameters:
131: +  ls - the SNES linesearch object
132: -  vf - the context for the monitor, in this case it is an ASCII PetscViewer and format

134:    Level: intermediate

136: .keywords: SNES, nonlinear, default, monitor, norm

138: .seealso: SNESGetLineSearch(), SNESMonitorSet(), SNESMonitorSolution()
139: @*/
140: PetscErrorCode  SNESLineSearchMonitorSolutionUpdate(SNESLineSearch ls,PetscViewerAndFormat *vf)
141: {
143:   PetscViewer    viewer = vf->viewer;
144:   Vec            Y,W,G;

147:   SNESLineSearchGetVecs(ls,NULL,NULL,&Y,&W,&G);
148:   PetscViewerPushFormat(viewer,vf->format);
149:   PetscViewerASCIIPrintf(viewer,"LineSearch attempted update to solution \n");
150:   VecView(Y,viewer);
151:   PetscViewerASCIIPrintf(viewer,"LineSearch attempted new solution \n");
152:   VecView(W,viewer);
153:   PetscViewerASCIIPrintf(viewer,"LineSearch attempted updated function value\n");
154:   VecView(G,viewer);
155:   PetscViewerPopFormat(viewer);
156:   return(0);
157: }

159: /*@
160:    SNESLineSearchCreate - Creates the line search context.

162:    Logically Collective on Comm

164:    Input Parameters:
165: .  comm - MPI communicator for the line search (typically from the associated SNES context).

167:    Output Parameters:
168: .  outlinesearch - the new linesearch context

170:    Level: developer

172:    Notes:
173:    The preferred calling sequence for users is to use SNESGetLineSearch() to acquire the SNESLineSearch instance
174:    already associated with the SNES.  This function is for developer use.

176: .keywords: LineSearch, create, context

178: .seealso: LineSearchDestroy(), SNESGetLineSearch()
179: @*/

181: PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch)
182: {
184:   SNESLineSearch linesearch;

188:   SNESInitializePackage();
189:   *outlinesearch = NULL;

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

193:   linesearch->vec_sol_new  = NULL;
194:   linesearch->vec_func_new = NULL;
195:   linesearch->vec_sol      = NULL;
196:   linesearch->vec_func     = NULL;
197:   linesearch->vec_update   = NULL;

199:   linesearch->lambda       = 1.0;
200:   linesearch->fnorm        = 1.0;
201:   linesearch->ynorm        = 1.0;
202:   linesearch->xnorm        = 1.0;
203:   linesearch->result       = SNES_LINESEARCH_SUCCEEDED;
204:   linesearch->norms        = PETSC_TRUE;
205:   linesearch->keeplambda   = PETSC_FALSE;
206:   linesearch->damping      = 1.0;
207:   linesearch->maxstep      = 1e8;
208:   linesearch->steptol      = 1e-12;
209:   linesearch->rtol         = 1e-8;
210:   linesearch->atol         = 1e-15;
211:   linesearch->ltol         = 1e-8;
212:   linesearch->precheckctx  = NULL;
213:   linesearch->postcheckctx = NULL;
214:   linesearch->max_its      = 1;
215:   linesearch->setupcalled  = PETSC_FALSE;
216:   *outlinesearch           = linesearch;
217:   return(0);
218: }

220: /*@
221:    SNESLineSearchSetUp - Prepares the line search for being applied by allocating
222:    any required vectors.

224:    Collective on SNESLineSearch

226:    Input Parameters:
227: .  linesearch - The LineSearch instance.

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

236:    Level: advanced

238: .keywords: SNESLineSearch, SetUp

240: .seealso: SNESGetLineSearch(), SNESLineSearchReset()
241: @*/

243: PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch)
244: {

248:   if (!((PetscObject)linesearch)->type_name) {
249:     SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);
250:   }
251:   if (!linesearch->setupcalled) {
252:     if (!linesearch->vec_sol_new) {
253:       VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new);
254:     }
255:     if (!linesearch->vec_func_new) {
256:       VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new);
257:     }
258:     if (linesearch->ops->setup) {
259:       (*linesearch->ops->setup)(linesearch);
260:     }
261:     if (!linesearch->ops->snesfunc) {SNESLineSearchSetFunction(linesearch,SNESComputeFunction);}
262:     linesearch->lambda      = linesearch->damping;
263:     linesearch->setupcalled = PETSC_TRUE;
264:   }
265:   return(0);
266: }


269: /*@
270:    SNESLineSearchReset - Undoes the SNESLineSearchSetUp() and deletes any Vecs or Mats allocated by the line search.

272:    Collective on SNESLineSearch

274:    Input Parameters:
275: .  linesearch - The LineSearch instance.

277:    Notes: Usually only called by SNESReset()

279:    Level: developer

281: .keywords: SNESLineSearch, Reset

283: .seealso: SNESGetLineSearch(), SNESLineSearchSetUp()
284: @*/

286: PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch)
287: {

291:   if (linesearch->ops->reset) (*linesearch->ops->reset)(linesearch);

293:   VecDestroy(&linesearch->vec_sol_new);
294:   VecDestroy(&linesearch->vec_func_new);

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

298:   linesearch->nwork       = 0;
299:   linesearch->setupcalled = PETSC_FALSE;
300:   return(0);
301: }

303: /*@C
304:    SNESLineSearchSetFunction - Sets the function evaluation used by the SNES line search

306:    Input Parameters:
307: .  linesearch - the SNESLineSearch context
308: +  func       - function evaluation routine

310:    Level: developer

312:    Notes: This is used internally by PETSc and not called by users

314: .keywords: get, linesearch, pre-check

316: .seealso: SNESGetLineSearch(), SNESSetFunction()
317: @*/
318: PetscErrorCode  SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES,Vec,Vec))
319: {
322:   linesearch->ops->snesfunc = func;
323:   return(0);
324: }

326: /*@C
327:    SNESLineSearchSetPreCheck - Sets a user function that is called after the initial search direction has been computed but
328:          before the line search routine has been applied. Allows the user to adjust the result of (usually a linear solve) that
329:          determined the search direction.

331:    Logically Collective on SNESLineSearch

333:    Input Parameters:
334: +  linesearch - the SNESLineSearch context
335: .  func - [optional] function evaluation routine, see SNESLineSearchPreCheck() for the calling sequence
336: -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)

338:    Level: intermediate

340: .keywords: set, linesearch, pre-check

342: .seealso: SNESGetLineSearch(), SNESLineSearchPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
343: @*/
344: PetscErrorCode  SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void *ctx)
345: {
348:   if (func) linesearch->ops->precheck = func;
349:   if (ctx) linesearch->precheckctx = ctx;
350:   return(0);
351: }

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

356:    Input Parameters:
357: .  linesearch - the SNESLineSearch context

359:    Output Parameters:
360: +  func       - [optional] function evaluation routine, see SNESLineSearchPreCheck() for calling sequence
361: -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)

363:    Level: intermediate

365: .keywords: get, linesearch, pre-check

367: .seealso: SNESGetLineSearch(), SNESGetLineSearch(), SNESLineSearchPreCheck(), SNESLineSearchGetPostCheck(), SNESLineSearchSetPreCheck(), SNESLineSearchSetPostCheck()
368: @*/
369: PetscErrorCode  SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void **ctx)
370: {
373:   if (func) *func = linesearch->ops->precheck;
374:   if (ctx) *ctx = linesearch->precheckctx;
375:   return(0);
376: }

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

382:    Logically Collective on SNESLineSearch

384:    Input Parameters:
385: +  linesearch - the SNESLineSearch context
386: .  func - [optional] function evaluation routine, see SNESLineSearchPostCheck()  for the calling sequence
387: -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)

389:    Level: intermediate

391: .keywords: set, linesearch, post-check

393: .seealso: SNESGetLineSearch(), SNESLineSearchPostCheck(), SNESLineSearchSetPreCheck(), SNESLineSearchGetPreCheck(), SNESLineSearchGetPostCheck()
394: @*/
395: PetscErrorCode  SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void *ctx)
396: {
399:   if (func) linesearch->ops->postcheck = func;
400:   if (ctx) linesearch->postcheckctx = ctx;
401:   return(0);
402: }

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

407:    Input Parameters:
408: .  linesearch - the SNESLineSearch context

410:    Output Parameters:
411: +  func - [optional] function evaluation routine, see for the calling sequence SNESLineSearchPostCheck()
412: -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)

414:    Level: intermediate

416: .keywords: get, linesearch, post-check

418: .seealso: SNESGetLineSearch(), SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchPostCheck(), SNESLineSearchSetPreCheck()
419: @*/
420: PetscErrorCode  SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void **ctx)
421: {
424:   if (func) *func = linesearch->ops->postcheck;
425:   if (ctx) *ctx = linesearch->postcheckctx;
426:   return(0);
427: }

429: /*@
430:    SNESLineSearchPreCheck - Prepares the line search for being applied.

432:    Logically Collective on SNESLineSearch

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

439:    Output Parameters:
440: .  changed - Indicator that the precheck routine has changed anything

442:    Level: developer

444: .keywords: SNESLineSearch, Create

446: .seealso: SNESGetLineSearch(), SNESLineSearchPostCheck(), SNESLineSearchSetPreCheck(), SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck()
447: @*/
448: PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed)
449: {

453:   *changed = PETSC_FALSE;
454:   if (linesearch->ops->precheck) {
455:     (*linesearch->ops->precheck)(linesearch, X, Y, changed, linesearch->precheckctx);
457:   }
458:   return(0);
459: }

461: /*@
462:    SNESLineSearchPostCheck - Prepares the line search for being applied.

464:    Logically Collective on SNESLineSearch

466:    Input Parameters:
467: +  linesearch - The linesearch context
468: .  X - The last solution
469: .  Y - The step direction
470: -  W - The updated solution, W = X + lambda*Y for some lambda

472:    Output Parameters:
473: +  changed_Y - Indicator if the direction Y has been changed.
474: -  changed_W - Indicator if the new candidate solution W has been changed.

476:    Level: developer

478: .keywords: SNESLineSearch, Create

480: .seealso: SNESGetLineSearch(), SNESLineSearchPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck(), SNESLineSearchSetPrecheck(), SNESLineSearchGetPrecheck()
481: @*/
482: PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch,Vec X,Vec Y,Vec W,PetscBool *changed_Y,PetscBool *changed_W)
483: {

487:   *changed_Y = PETSC_FALSE;
488:   *changed_W = PETSC_FALSE;
489:   if (linesearch->ops->postcheck) {
490:     (*linesearch->ops->postcheck)(linesearch,X,Y,W,changed_Y,changed_W,linesearch->postcheckctx);
493:   }
494:   return(0);
495: }

497: /*@C
498:    SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration

500:    Logically Collective on SNESLineSearch

502:    Input Arguments:
503: +  linesearch - linesearch context
504: .  X - base state for this step
505: .  Y - initial correction
506: -  ctx - context for this function

508:    Output Arguments:
509: +  Y - correction, possibly modified
510: -  changed - flag indicating that Y was modified

512:    Options Database Key:
513: +  -snes_linesearch_precheck_picard - activate this routine
514: -  -snes_linesearch_precheck_picard_angle - angle

516:    Level: advanced

518:    Notes:
519:    This function should be passed to SNESLineSearchSetPreCheck()

521:    The justification for this method involves the linear convergence of a Picard iteration
522:    so the Picard linearization should be provided in place of the "Jacobian". This correction
523:    is generally not useful when using a Newton linearization.

525:    Reference:
526:    Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology.

528: .seealso: SNESGetLineSearch(), SNESLineSearchSetPreCheck()
529: @*/
530: PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed,void *ctx)
531: {
533:   PetscReal      angle = *(PetscReal*)linesearch->precheckctx;
534:   Vec            Ylast;
535:   PetscScalar    dot;
536:   PetscInt       iter;
537:   PetscReal      ynorm,ylastnorm,theta,angle_radians;
538:   SNES           snes;

541:   SNESLineSearchGetSNES(linesearch, &snes);
542:   PetscObjectQuery((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject*)&Ylast);
543:   if (!Ylast) {
544:     VecDuplicate(Y,&Ylast);
545:     PetscObjectCompose((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject)Ylast);
546:     PetscObjectDereference((PetscObject)Ylast);
547:   }
548:   SNESGetIterationNumber(snes,&iter);
549:   if (iter < 2) {
550:     VecCopy(Y,Ylast);
551:     *changed = PETSC_FALSE;
552:     return(0);
553:   }

555:   VecDot(Y,Ylast,&dot);
556:   VecNorm(Y,NORM_2,&ynorm);
557:   VecNorm(Ylast,NORM_2,&ylastnorm);
558:   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
559:   theta         = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0));
560:   angle_radians = angle * PETSC_PI / 180.;
561:   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
562:     /* Modify the step Y */
563:     PetscReal alpha,ydiffnorm;
564:     VecAXPY(Ylast,-1.0,Y);
565:     VecNorm(Ylast,NORM_2,&ydiffnorm);
566:     alpha = ylastnorm / ydiffnorm;
567:     VecCopy(Y,Ylast);
568:     VecScale(Y,alpha);
569:     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);
570:   } else {
571:     PetscInfo2(snes,"Angle %14.12e degrees exceeds threshold %14.12e, no correction applied\n",(double)(theta*180./PETSC_PI),(double)angle);
572:     VecCopy(Y,Ylast);
573:     *changed = PETSC_FALSE;
574:   }
575:   return(0);
576: }

578: /*@
579:    SNESLineSearchApply - Computes the line-search update.

581:    Collective on SNESLineSearch

583:    Input Parameters:
584: +  linesearch - The linesearch context
585: .  X - The current solution
586: .  F - The current function
587: .  fnorm - The current norm
588: -  Y - The search direction

590:    Output Parameters:
591: +  X - The new solution
592: .  F - The new function
593: -  fnorm - The new function norm

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

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

611:    Level: Intermediate

613: .keywords: SNESLineSearch, Create

615: .seealso: SNESGetLineSearch(), SNESLineSearchCreate(), SNESLineSearchPreCheck(), SNESLineSearchPostCheck(), SNESSolve(), SNESComputeFunction(), SNESLineSearchSetComputeNorms(),
616:           SNESLineSearchType, SNESLineSearchSetType()
617: @*/
618: PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y)
619: {


628:   linesearch->result = SNES_LINESEARCH_SUCCEEDED;

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

634:   SNESLineSearchSetUp(linesearch);

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

638:   if (fnorm) linesearch->fnorm = *fnorm;
639:   else {
640:     VecNorm(F, NORM_2, &linesearch->fnorm);
641:   }

643:   PetscLogEventBegin(SNESLINESEARCH_Apply,linesearch,X,F,Y);

645:   (*linesearch->ops->apply)(linesearch);

647:   PetscLogEventEnd(SNESLINESEARCH_Apply,linesearch,X,F,Y);

649:   if (fnorm) *fnorm = linesearch->fnorm;
650:   return(0);
651: }

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

656:    Collective on SNESLineSearch

658:    Input Parameters:
659: .  linesearch - The linesearch context

661:    Level: developer

663: .keywords: SNESLineSearch, Destroy

665: .seealso: SNESGetLineSearch(), SNESLineSearchCreate(), SNESLineSearchReset(), SNESDestroy()
666: @*/
667: PetscErrorCode SNESLineSearchDestroy(SNESLineSearch * linesearch)
668: {

672:   if (!*linesearch) return(0);
674:   if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = 0; return(0);}
675:   PetscObjectSAWsViewOff((PetscObject)*linesearch);
676:   SNESLineSearchReset(*linesearch);
677:   if ((*linesearch)->ops->destroy) (*linesearch)->ops->destroy(*linesearch);
678:   PetscViewerDestroy(&(*linesearch)->monitor);
679:   SNESLineSearchMonitorCancel((*linesearch));
680:   PetscHeaderDestroy(linesearch);
681:   return(0);
682: }

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

687:    Input Parameters:
688: +  linesearch - the linesearch object
689: -  viewer - an ASCII PetscViewer or NULL to turn off monitor

691:    Logically Collective on SNESLineSearch

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

696:    Level: intermediate

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

702: .seealso: SNESGetLineSearch(), SNESLineSearchGetDefaultMonitor(), PetscViewer, SNESLineSearchSetMonitor()
703: @*/
704: PetscErrorCode  SNESLineSearchSetDefaultMonitor(SNESLineSearch linesearch, PetscViewer viewer)
705: {

709:   if (viewer) {PetscObjectReference((PetscObject)viewer);}
710:   PetscViewerDestroy(&linesearch->monitor);
711:   linesearch->monitor = viewer;
712:   return(0);
713: }

715: /*@
716:    SNESLineSearchGetDefaultMonitor - Gets the PetscViewer instance for the line search monitor.

718:    Input Parameter:
719: .  linesearch - linesearch context

721:    Output Parameter:
722: .  monitor - monitor context

724:    Logically Collective on SNES

726:    Options Database Keys:
727: .   -snes_linesearch_monitor - enables the monitor

729:    Level: intermediate

731: .seealso: SNESGetLineSearch(), SNESLineSearchSetDefaultMonitor(), PetscViewer
732: @*/
733: PetscErrorCode  SNESLineSearchGetDefaultMonitor(SNESLineSearch linesearch, PetscViewer *monitor)
734: {
737:   if (monitor) {
739:     *monitor = linesearch->monitor;
740:   }
741:   return(0);
742: }

744: /*@C
745:    SNESLineSearchMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user

747:    Collective on SNESLineSearch

749:    Input Parameters:
750: +  ls - LineSearch object you wish to monitor
751: .  name - the monitor type one is seeking
752: .  help - message indicating what monitoring is done
753: .  manual - manual page for the monitor
754: .  monitor - the monitor function
755: -  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

757:    Level: developer

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

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

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

791:    Input Parameters:
792: .  linesearch - linesearch context

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

811:    Logically Collective on SNESLineSearch

813:    Level: intermediate

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

827:   SNESLineSearchRegisterAll();

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

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

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

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

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

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

864:       PetscOptionsReal("-snes_linesearch_precheck_picard_angle","Maximum angle at which to activate the correction",
865:                               "none",linesearch->precheck_picard_angle,&linesearch->precheck_picard_angle,NULL);
866:       SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&linesearch->precheck_picard_angle);
867:     } else {
868:       SNESLineSearchSetPreCheck(linesearch,NULL,NULL);
869:     }
870:   }
871:   PetscOptionsInt("-snes_linesearch_order","Order of approximation used in the line search","SNESLineSearchSetOrder",linesearch->order,&linesearch->order,NULL);
872:   PetscOptionsBool("-snes_linesearch_norms","Compute final norms in line search","SNESLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,NULL);

874:   if (linesearch->ops->setfromoptions) {
875:     (*linesearch->ops->setfromoptions)(PetscOptionsObject,linesearch);
876:   }

878:   PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)linesearch);
879:   PetscOptionsEnd();
880:   return(0);
881: }

883: /*@
884:    SNESLineSearchView - Prints useful information about the line search

886:    Input Parameters:
887: .  linesearch - linesearch context

889:    Logically Collective on SNESLineSearch

891:    Level: intermediate

893: .seealso: SNESLineSearchCreate()
894: @*/
895: PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer)
896: {
898:   PetscBool      iascii;

902:   if (!viewer) {
903:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch),&viewer);
904:   }

908:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
909:   if (iascii) {
910:     PetscObjectPrintClassNamePrefixType((PetscObject)linesearch,viewer);
911:     if (linesearch->ops->view) {
912:       PetscViewerASCIIPushTab(viewer);
913:       (*linesearch->ops->view)(linesearch,viewer);
914:       PetscViewerASCIIPopTab(viewer);
915:     }
916:     PetscViewerASCIIPrintf(viewer,"  maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep,(double)linesearch->steptol);
917:     PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol,(double)linesearch->atol,(double)linesearch->ltol);
918:     PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D\n", linesearch->max_its);
919:     if (linesearch->ops->precheck) {
920:       if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) {
921:         PetscViewerASCIIPrintf(viewer,"  using precheck step to speed up Picard convergence\n", linesearch->max_its);
922:       } else {
923:         PetscViewerASCIIPrintf(viewer,"  using user-defined precheck step\n", linesearch->max_its);
924:       }
925:     }
926:     if (linesearch->ops->postcheck) {
927:       PetscViewerASCIIPrintf(viewer,"  using user-defined postcheck step\n", linesearch->max_its);
928:     }
929:   }
930:   return(0);
931: }

933: /*@C
934:    SNESLineSearchSetType - Sets the linesearch type

936:    Logically Collective on SNESLineSearch

938:    Input Parameters:
939: +  linesearch - linesearch context
940: -  type - The type of line search to be used

942:    Available Types:
943: +  SNESLINESEARCHBASIC - Simple damping line search, defaults to using the full Newton step
944: .  SNESLINESEARCHBT - Backtracking line search over the L2 norm of the function
945: .  SNESLINESEARCHL2 - Secant line search over the L2 norm of the function
946: .  SNESLINESEARCHCP - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x)
947: .  SNESLINESEARCHNLEQERR - Affine-covariant error-oriented linesearch
948: -  SNESLINESEARCHSHELL - User provided SNESLineSearch implementation

950:    Options Database:
951: .  -snes_linesearch_type <type> - basic, bt, l2, cp, nleqerr, shell

953:    Level: intermediate

955: .seealso: SNESLineSearchCreate(), SNESLineSearchType, SNESLineSearchSetFromOptions()
956: @*/
957: PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type)
958: {
959:   PetscErrorCode ierr,(*r)(SNESLineSearch);
960:   PetscBool      match;


966:   PetscObjectTypeCompare((PetscObject)linesearch,type,&match);
967:   if (match) return(0);

969:   PetscFunctionListFind(SNESLineSearchList,type,&r);
970:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type);
971:   /* Destroy the previous private linesearch context */
972:   if (linesearch->ops->destroy) {
973:     (*(linesearch)->ops->destroy)(linesearch);

975:     linesearch->ops->destroy = NULL;
976:   }
977:   /* Reinitialize function pointers in SNESLineSearchOps structure */
978:   linesearch->ops->apply          = 0;
979:   linesearch->ops->view           = 0;
980:   linesearch->ops->setfromoptions = 0;
981:   linesearch->ops->destroy        = 0;

983:   PetscObjectChangeTypeName((PetscObject)linesearch,type);
984:   (*r)(linesearch);
985:   return(0);
986: }

988: /*@
989:    SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation.

991:    Input Parameters:
992: +  linesearch - linesearch context
993: -  snes - The snes instance

995:    Level: developer

997:    Notes:
998:    This happens automatically when the line search is obtained/created with
999:    SNESGetLineSearch().  This routine is therefore mainly called within SNES
1000:    implementations.

1002:    Level: developer

1004: .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
1005: @*/
1006: PetscErrorCode  SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes)
1007: {
1011:   linesearch->snes = snes;
1012:   return(0);
1013: }

1015: /*@
1016:    SNESLineSearchGetSNES - Gets the SNES instance associated with the line search.
1017:    Having an associated SNES is necessary because most line search implementations must be able to
1018:    evaluate the function using SNESComputeFunction() for the associated SNES.  This routine
1019:    is used in the line search implementations when one must get this associated SNES instance.

1021:    Input Parameters:
1022: .  linesearch - linesearch context

1024:    Output Parameters:
1025: .  snes - The snes instance

1027:    Level: developer

1029: .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
1030: @*/
1031: PetscErrorCode  SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes)
1032: {
1036:   *snes = linesearch->snes;
1037:   return(0);
1038: }

1040: /*@
1041:    SNESLineSearchGetLambda - Gets the last linesearch steplength discovered.

1043:    Input Parameters:
1044: .  linesearch - linesearch context

1046:    Output Parameters:
1047: .  lambda - The last steplength computed during SNESLineSearchApply()

1049:    Level: advanced

1051:    Notes:
1052:    This is useful in methods where the solver is ill-scaled and
1053:    requires some adaptive notion of the difference in scale between the
1054:    solution and the function.  For instance, SNESQN may be scaled by the
1055:    line search lambda using the argument -snes_qn_scaling ls.

1057: .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetDamping(), SNESLineSearchApply()
1058: @*/
1059: PetscErrorCode  SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda)
1060: {
1064:   *lambda = linesearch->lambda;
1065:   return(0);
1066: }

1068: /*@
1069:    SNESLineSearchSetLambda - Sets the linesearch steplength.

1071:    Input Parameters:
1072: +  linesearch - linesearch context
1073: -  lambda - The last steplength.

1075:    Notes:
1076:    This routine is typically used within implementations of SNESLineSearchApply()
1077:    to set the final steplength.  This routine (and SNESLineSearchGetLambda()) were
1078:    added in order to facilitate Quasi-Newton methods that use the previous steplength
1079:    as an inner scaling parameter.

1081:    Level: advanced

1083: .seealso: SNESLineSearchGetLambda()
1084: @*/
1085: PetscErrorCode  SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
1086: {
1089:   linesearch->lambda = lambda;
1090:   return(0);
1091: }

1093: /*@
1094:    SNESLineSearchGetTolerances - Gets the tolerances for the linesearch.  These include
1095:    tolerances for the relative and absolute change in the function norm, the change
1096:    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
1097:    and the maximum number of iterations the line search procedure may take.

1099:    Input Parameters:
1100: .  linesearch - linesearch context

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

1110:    Level: intermediate

1112:    Notes:
1113:    Different line searches may implement these parameters slightly differently as
1114:    the type requires.

1116: .seealso: SNESLineSearchSetTolerances()
1117: @*/
1118: PetscErrorCode  SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
1119: {
1122:   if (steptol) {
1124:     *steptol = linesearch->steptol;
1125:   }
1126:   if (maxstep) {
1128:     *maxstep = linesearch->maxstep;
1129:   }
1130:   if (rtol) {
1132:     *rtol = linesearch->rtol;
1133:   }
1134:   if (atol) {
1136:     *atol = linesearch->atol;
1137:   }
1138:   if (ltol) {
1140:     *ltol = linesearch->ltol;
1141:   }
1142:   if (max_its) {
1144:     *max_its = linesearch->max_its;
1145:   }
1146:   return(0);
1147: }

1149: /*@
1150:    SNESLineSearchSetTolerances -  Gets the tolerances for the linesearch.  These include
1151:    tolerances for the relative and absolute change in the function norm, the change
1152:    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
1153:    and the maximum number of iterations the line search procedure may take.

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

1164:    Notes:
1165:    The user may choose to not set any of the tolerances using PETSC_DEFAULT in
1166:    place of an argument.

1168:    Level: intermediate

1170: .seealso: SNESLineSearchGetTolerances()
1171: @*/
1172: PetscErrorCode  SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its)
1173: {

1183:   if (steptol!= PETSC_DEFAULT) {
1184:     if (steptol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Minimum step length %14.12e must be non-negative",(double)steptol);
1185:     linesearch->steptol = steptol;
1186:   }

1188:   if (maxstep!= PETSC_DEFAULT) {
1189:     if (maxstep < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum step length %14.12e must be non-negative",(double)maxstep);
1190:     linesearch->maxstep = maxstep;
1191:   }

1193:   if (rtol != PETSC_DEFAULT) {
1194:     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);
1195:     linesearch->rtol = rtol;
1196:   }

1198:   if (atol != PETSC_DEFAULT) {
1199:     if (atol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %14.12e must be non-negative",(double)atol);
1200:     linesearch->atol = atol;
1201:   }

1203:   if (ltol != PETSC_DEFAULT) {
1204:     if (ltol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Labmda tolerance %14.12e must be non-negative",(double)ltol);
1205:     linesearch->ltol = ltol;
1206:   }

1208:   if (max_its != PETSC_DEFAULT) {
1209:     if (max_its < 0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",max_its);
1210:     linesearch->max_its = max_its;
1211:   }
1212:   return(0);
1213: }

1215: /*@
1216:    SNESLineSearchGetDamping - Gets the line search damping parameter.

1218:    Input Parameters:
1219: .  linesearch - linesearch context

1221:    Output Parameters:
1222: .  damping - The damping parameter

1224:    Level: advanced

1226: .seealso: SNESLineSearchGetStepTolerance(), SNESQN
1227: @*/

1229: PetscErrorCode  SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping)
1230: {
1234:   *damping = linesearch->damping;
1235:   return(0);
1236: }

1238: /*@
1239:    SNESLineSearchSetDamping - Sets the line search damping paramter.

1241:    Input Parameters:
1242: +  linesearch - linesearch context
1243: -  damping - The damping parameter

1245:    Options Database:
1246: .   -snes_linesearch_damping
1247:    Level: intermediate

1249:    Notes:
1250:    The basic line search merely takes the update step scaled by the damping parameter.
1251:    The use of the damping parameter in the l2 and cp line searches is much more subtle;
1252:    it is used as a starting point in calculating the secant step. However, the eventual
1253:    step may be of greater length than the damping parameter.  In the bt line search it is
1254:    used as the maximum possible step length, as the bt line search only backtracks.

1256: .seealso: SNESLineSearchGetDamping()
1257: @*/
1258: PetscErrorCode  SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping)
1259: {
1262:   linesearch->damping = damping;
1263:   return(0);
1264: }

1266: /*@
1267:    SNESLineSearchGetOrder - Gets the line search approximation order.

1269:    Input Parameters:
1270: .  linesearch - linesearch context

1272:    Output Parameters:
1273: .  order - The order

1275:    Possible Values for order:
1276: +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
1277: .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
1278: -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order

1280:    Level: intermediate

1282: .seealso: SNESLineSearchSetOrder()
1283: @*/

1285: PetscErrorCode  SNESLineSearchGetOrder(SNESLineSearch linesearch,PetscInt *order)
1286: {
1290:   *order = linesearch->order;
1291:   return(0);
1292: }

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

1297:    Input Parameters:
1298: .  linesearch - linesearch context
1299: .  order - The damping parameter

1301:    Level: intermediate

1303:    Possible Values for order:
1304: +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
1305: .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
1306: -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order

1308:    Notes:
1309:    Variable orders are supported by the following line searches:
1310: +  bt - cubic and quadratic
1311: -  cp - linear and quadratic

1313: .seealso: SNESLineSearchGetOrder(), SNESLineSearchSetDamping()
1314: @*/
1315: PetscErrorCode  SNESLineSearchSetOrder(SNESLineSearch linesearch,PetscInt order)
1316: {
1319:   linesearch->order = order;
1320:   return(0);
1321: }

1323: /*@
1324:    SNESLineSearchGetNorms - Gets the norms for for X, Y, and F.

1326:    Input Parameters:
1327: .  linesearch - linesearch context

1329:    Output Parameters:
1330: +  xnorm - The norm of the current solution
1331: .  fnorm - The norm of the current function
1332: -  ynorm - The norm of the current update

1334:    Notes:
1335:    This function is mainly called from SNES implementations.

1337:    Level: developer

1339: .seealso: SNESLineSearchSetNorms() SNESLineSearchGetVecs()
1340: @*/
1341: PetscErrorCode  SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm)
1342: {
1345:   if (xnorm) *xnorm = linesearch->xnorm;
1346:   if (fnorm) *fnorm = linesearch->fnorm;
1347:   if (ynorm) *ynorm = linesearch->ynorm;
1348:   return(0);
1349: }

1351: /*@
1352:    SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F.

1354:    Input Parameters:
1355: +  linesearch - linesearch context
1356: .  xnorm - The norm of the current solution
1357: .  fnorm - The norm of the current function
1358: -  ynorm - The norm of the current update

1360:    Level: advanced

1362: .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1363: @*/
1364: PetscErrorCode  SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
1365: {
1368:   linesearch->xnorm = xnorm;
1369:   linesearch->fnorm = fnorm;
1370:   linesearch->ynorm = ynorm;
1371:   return(0);
1372: }

1374: /*@
1375:    SNESLineSearchComputeNorms - Computes the norms of X, F, and Y.

1377:    Input Parameters:
1378: .  linesearch - linesearch context

1380:    Options Database Keys:
1381: .   -snes_linesearch_norms - turn norm computation on or off

1383:    Level: intermediate

1385: .seealso: SNESLineSearchGetNorms, SNESLineSearchSetNorms(), SNESLineSearchSetComputeNorms()
1386: @*/
1387: PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
1388: {
1390:   SNES           snes;

1393:   if (linesearch->norms) {
1394:     if (linesearch->ops->vinorm) {
1395:       SNESLineSearchGetSNES(linesearch, &snes);
1396:       VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);
1397:       VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);
1398:       (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);
1399:     } else {
1400:       VecNormBegin(linesearch->vec_func,   NORM_2, &linesearch->fnorm);
1401:       VecNormBegin(linesearch->vec_sol,    NORM_2, &linesearch->xnorm);
1402:       VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);
1403:       VecNormEnd(linesearch->vec_func,     NORM_2, &linesearch->fnorm);
1404:       VecNormEnd(linesearch->vec_sol,      NORM_2, &linesearch->xnorm);
1405:       VecNormEnd(linesearch->vec_update,   NORM_2, &linesearch->ynorm);
1406:     }
1407:   }
1408:   return(0);
1409: }

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

1414:    Input Parameters:
1415: +  linesearch  - linesearch context
1416: -  flg  - indicates whether or not to compute norms

1418:    Options Database Keys:
1419: .   -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic linesearch

1421:    Notes:
1422:    This is most relevant to the SNESLINESEARCHBASIC line search type since most line searches have a stopping criteria involving the norm.

1424:    Level: intermediate

1426: .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetNorms(), SNESLineSearchComputeNorms(), SNESLINESEARCHBASIC
1427: @*/
1428: PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg)
1429: {
1431:   linesearch->norms = flg;
1432:   return(0);
1433: }

1435: /*@
1436:    SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context

1438:    Input Parameters:
1439: .  linesearch - linesearch context

1441:    Output Parameters:
1442: +  X - Solution vector
1443: .  F - Function vector
1444: .  Y - Search direction vector
1445: .  W - Solution work vector
1446: -  G - Function work vector

1448:    Notes:
1449:    At the beginning of a line search application, X should contain a
1450:    solution and the vector F the function computed at X.  At the end of the
1451:    line search application, X should contain the new solution, and F the
1452:    function evaluated at the new solution.

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

1456:    Level: advanced

1458: .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1459: @*/
1460: PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G)
1461: {
1464:   if (X) {
1466:     *X = linesearch->vec_sol;
1467:   }
1468:   if (F) {
1470:     *F = linesearch->vec_func;
1471:   }
1472:   if (Y) {
1474:     *Y = linesearch->vec_update;
1475:   }
1476:   if (W) {
1478:     *W = linesearch->vec_sol_new;
1479:   }
1480:   if (G) {
1482:     *G = linesearch->vec_func_new;
1483:   }
1484:   return(0);
1485: }

1487: /*@
1488:    SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context

1490:    Input Parameters:
1491: +  linesearch - linesearch context
1492: .  X - Solution vector
1493: .  F - Function vector
1494: .  Y - Search direction vector
1495: .  W - Solution work vector
1496: -  G - Function work vector

1498:    Level: advanced

1500: .seealso: SNESLineSearchSetNorms(), SNESLineSearchGetVecs()
1501: @*/
1502: PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G)
1503: {
1506:   if (X) {
1508:     linesearch->vec_sol = X;
1509:   }
1510:   if (F) {
1512:     linesearch->vec_func = F;
1513:   }
1514:   if (Y) {
1516:     linesearch->vec_update = Y;
1517:   }
1518:   if (W) {
1520:     linesearch->vec_sol_new = W;
1521:   }
1522:   if (G) {
1524:     linesearch->vec_func_new = G;
1525:   }
1526:   return(0);
1527: }

1529: /*@C
1530:    SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1531:    SNES options in the database.

1533:    Logically Collective on SNESLineSearch

1535:    Input Parameters:
1536: +  snes - the SNES context
1537: -  prefix - the prefix to prepend to all option names

1539:    Notes:
1540:    A hyphen (-) must NOT be given at the beginning of the prefix name.
1541:    The first character of all runtime options is AUTOMATICALLY the hyphen.

1543:    Level: advanced

1545: .keywords: SNESLineSearch, append, options, prefix, database

1547: .seealso: SNESGetOptionsPrefix()
1548: @*/
1549: PetscErrorCode  SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[])
1550: {

1555:   PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);
1556:   return(0);
1557: }

1559: /*@C
1560:    SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all
1561:    SNESLineSearch options in the database.

1563:    Not Collective

1565:    Input Parameter:
1566: .  linesearch - the SNESLineSearch context

1568:    Output Parameter:
1569: .  prefix - pointer to the prefix string used

1571:    Notes:
1572:    On the fortran side, the user should pass in a string 'prefix' of
1573:    sufficient length to hold the prefix.

1575:    Level: advanced

1577: .keywords: SNESLineSearch, get, options, prefix, database

1579: .seealso: SNESAppendOptionsPrefix()
1580: @*/
1581: PetscErrorCode  SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[])
1582: {

1587:   PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);
1588:   return(0);
1589: }

1591: /*@C
1592:    SNESLineSearchSetWorkVecs - Gets work vectors for the line search.

1594:    Input Parameter:
1595: +  linesearch - the SNESLineSearch context
1596: -  nwork - the number of work vectors

1598:    Level: developer

1600:    Developers Note: This is PETSC_EXTERN because it may be used by user written plugin SNES implementations

1602: .keywords: SNESLineSearch, work, vector

1604: .seealso: SNESSetWorkVecs()
1605: @*/
1606: PetscErrorCode  SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork)
1607: {

1611:   if (linesearch->vec_sol) {
1612:     VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);
1613:   } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
1614:   return(0);
1615: }

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

1620:    Input Parameters:
1621: .  linesearch - linesearch context

1623:    Output Parameters:
1624: .  result - The success or failure status

1626:    Notes:
1627:    This is typically called after SNESLineSearchApply() in order to determine if the line-search failed
1628:    (and set the SNES convergence accordingly).

1630:    Level: intermediate

1632: .seealso: SNESLineSearchSetReason(), SNESLineSearchReason
1633: @*/
1634: PetscErrorCode  SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result)
1635: {
1639:   *result = linesearch->result;
1640:   return(0);
1641: }

1643: /*@
1644:    SNESLineSearchSetReason - Sets the success/failure status of the last line search application

1646:    Input Parameters:
1647: +  linesearch - linesearch context
1648: -  result - The success or failure status

1650:    Notes:
1651:    This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set
1652:    the success or failure of the line search method.

1654:    Level: developer

1656: .seealso: SNESLineSearchGetSResult()
1657: @*/
1658: PetscErrorCode  SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result)
1659: {
1662:   linesearch->result = result;
1663:   return(0);
1664: }

1666: /*@C
1667:    SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.

1669:    Input Parameters:
1670: +  snes - nonlinear context obtained from SNESCreate()
1671: .  projectfunc - function for projecting the function to the bounds
1672: -  normfunc - function for computing the norm of an active set

1674:    Logically Collective on SNES

1676:    Calling sequence of projectfunc:
1677: .vb
1678:    projectfunc (SNES snes, Vec X)
1679: .ve

1681:     Input parameters for projectfunc:
1682: +   snes - nonlinear context
1683: -   X - current solution

1685:     Output parameters for projectfunc:
1686: .   X - Projected solution

1688:    Calling sequence of normfunc:
1689: .vb
1690:    projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm)
1691: .ve

1693:     Input parameters for normfunc:
1694: +   snes - nonlinear context
1695: .   X - current solution
1696: -   F - current residual

1698:     Output parameters for normfunc:
1699: .   fnorm - VI-specific norm of the function

1701:     Notes:
1702:     The VI solvers require projection of the solution to the feasible set.  projectfunc should implement this.

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

1707:     Level: developer

1709: .keywords: SNES, line search, VI, nonlinear, set, line search

1711: .seealso: SNESLineSearchGetVIFunctions(), SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck()
1712: @*/
1713: extern PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc)
1714: {
1717:   if (projectfunc) linesearch->ops->viproject = projectfunc;
1718:   if (normfunc) linesearch->ops->vinorm = normfunc;
1719:   return(0);
1720: }

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

1725:    Input Parameters:
1726: .  linesearch - the line search context, obtain with SNESGetLineSearch()

1728:    Output Parameters:
1729: +  projectfunc - function for projecting the function to the bounds
1730: -  normfunc - function for computing the norm of an active set

1732:    Logically Collective on SNES

1734:     Level: developer

1736: .keywords: SNES, line search, VI, nonlinear, get, line search

1738: .seealso: SNESLineSearchSetVIFunctions(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
1739: @*/
1740: extern PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc)
1741: {
1743:   if (projectfunc) *projectfunc = linesearch->ops->viproject;
1744:   if (normfunc) *normfunc = linesearch->ops->vinorm;
1745:   return(0);
1746: }

1748: /*@C
1749:   SNESLineSearchRegister - See SNESLineSearchRegister()

1751:   Level: advanced
1752: @*/
1753: PetscErrorCode  SNESLineSearchRegister(const char sname[],PetscErrorCode (*function)(SNESLineSearch))
1754: {

1758:   PetscFunctionListAdd(&SNESLineSearchList,sname,function);
1759:   return(0);
1760: }